1be1d678aSKris Buschelman #define PETSCMAT_DLL 2be1d678aSKris Buschelman 3ed3cc1f0SBarry Smith /* 4ed3cc1f0SBarry Smith Basic functions for basic parallel dense matrices. 5ed3cc1f0SBarry Smith */ 6ed3cc1f0SBarry Smith 770f55243SBarry Smith #include "src/mat/impls/dense/mpi/mpidense.h" 88965ea79SLois Curfman McInnes 9ba8c8a56SBarry Smith #undef __FUNCT__ 10ba8c8a56SBarry Smith #define __FUNCT__ "MatGetRow_MPIDense" 11ba8c8a56SBarry Smith PetscErrorCode MatGetRow_MPIDense(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v) 12ba8c8a56SBarry Smith { 13ba8c8a56SBarry Smith Mat_MPIDense *mat = (Mat_MPIDense*)A->data; 14ba8c8a56SBarry Smith PetscErrorCode ierr; 15ba8c8a56SBarry Smith PetscInt lrow,rstart = mat->rstart,rend = mat->rend; 16ba8c8a56SBarry Smith 17ba8c8a56SBarry Smith PetscFunctionBegin; 18ba8c8a56SBarry Smith if (row < rstart || row >= rend) SETERRQ(PETSC_ERR_SUP,"only local rows") 19ba8c8a56SBarry Smith lrow = row - rstart; 20ba8c8a56SBarry Smith ierr = MatGetRow(mat->A,lrow,nz,(const PetscInt **)idx,(const PetscScalar **)v);CHKERRQ(ierr); 21ba8c8a56SBarry Smith PetscFunctionReturn(0); 22ba8c8a56SBarry Smith } 23ba8c8a56SBarry Smith 24ba8c8a56SBarry Smith #undef __FUNCT__ 25ba8c8a56SBarry Smith #define __FUNCT__ "MatRestoreRow_MPIDense" 26ba8c8a56SBarry Smith PetscErrorCode MatRestoreRow_MPIDense(Mat mat,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v) 27ba8c8a56SBarry Smith { 28ba8c8a56SBarry Smith PetscErrorCode ierr; 29ba8c8a56SBarry Smith 30ba8c8a56SBarry Smith PetscFunctionBegin; 31ba8c8a56SBarry Smith if (idx) {ierr = PetscFree(*idx);CHKERRQ(ierr);} 32ba8c8a56SBarry Smith if (v) {ierr = PetscFree(*v);CHKERRQ(ierr);} 33ba8c8a56SBarry Smith PetscFunctionReturn(0); 34ba8c8a56SBarry Smith } 35ba8c8a56SBarry Smith 360de54da6SSatish Balay EXTERN_C_BEGIN 374a2ae208SSatish Balay #undef __FUNCT__ 384a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonalBlock_MPIDense" 39be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatGetDiagonalBlock_MPIDense(Mat A,PetscTruth *iscopy,MatReuse reuse,Mat *B) 400de54da6SSatish Balay { 410de54da6SSatish Balay Mat_MPIDense *mdn = (Mat_MPIDense*)A->data; 426849ba73SBarry Smith PetscErrorCode ierr; 4313f74950SBarry Smith PetscInt m = A->m,rstart = mdn->rstart; 4487828ca2SBarry Smith PetscScalar *array; 450de54da6SSatish Balay MPI_Comm comm; 460de54da6SSatish Balay 470de54da6SSatish Balay PetscFunctionBegin; 48273d9f13SBarry Smith if (A->M != A->N) SETERRQ(PETSC_ERR_SUP,"Only square matrices supported."); 490de54da6SSatish Balay 500de54da6SSatish Balay /* The reuse aspect is not implemented efficiently */ 510de54da6SSatish Balay if (reuse) { ierr = MatDestroy(*B);CHKERRQ(ierr);} 520de54da6SSatish Balay 530de54da6SSatish Balay ierr = PetscObjectGetComm((PetscObject)(mdn->A),&comm);CHKERRQ(ierr); 540de54da6SSatish Balay ierr = MatGetArray(mdn->A,&array);CHKERRQ(ierr); 55f69a0ea3SMatthew Knepley ierr = MatCreate(comm,B);CHKERRQ(ierr); 56f69a0ea3SMatthew Knepley ierr = MatSetSizes(*B,m,m,m,m);CHKERRQ(ierr); 575c5985e7SKris Buschelman ierr = MatSetType(*B,mdn->A->type_name);CHKERRQ(ierr); 585c5985e7SKris Buschelman ierr = MatSeqDenseSetPreallocation(*B,array+m*rstart);CHKERRQ(ierr); 590de54da6SSatish Balay ierr = MatRestoreArray(mdn->A,&array);CHKERRQ(ierr); 600de54da6SSatish Balay ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 610de54da6SSatish Balay ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 620de54da6SSatish Balay 630de54da6SSatish Balay *iscopy = PETSC_TRUE; 640de54da6SSatish Balay PetscFunctionReturn(0); 650de54da6SSatish Balay } 660de54da6SSatish Balay EXTERN_C_END 670de54da6SSatish Balay 684a2ae208SSatish Balay #undef __FUNCT__ 694a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_MPIDense" 7013f74950SBarry Smith PetscErrorCode MatSetValues_MPIDense(Mat mat,PetscInt m,const PetscInt idxm[],PetscInt n,const PetscInt idxn[],const PetscScalar v[],InsertMode addv) 718965ea79SLois Curfman McInnes { 7239b7565bSBarry Smith Mat_MPIDense *A = (Mat_MPIDense*)mat->data; 73dfbe8321SBarry Smith PetscErrorCode ierr; 7413f74950SBarry Smith PetscInt i,j,rstart = A->rstart,rend = A->rend,row; 75273d9f13SBarry Smith PetscTruth roworiented = A->roworiented; 768965ea79SLois Curfman McInnes 773a40ed3dSBarry Smith PetscFunctionBegin; 788965ea79SLois Curfman McInnes for (i=0; i<m; i++) { 795ef9f2a5SBarry Smith if (idxm[i] < 0) continue; 80273d9f13SBarry Smith if (idxm[i] >= mat->M) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large"); 818965ea79SLois Curfman McInnes if (idxm[i] >= rstart && idxm[i] < rend) { 828965ea79SLois Curfman McInnes row = idxm[i] - rstart; 8339b7565bSBarry Smith if (roworiented) { 8439b7565bSBarry Smith ierr = MatSetValues(A->A,1,&row,n,idxn,v+i*n,addv);CHKERRQ(ierr); 853a40ed3dSBarry Smith } else { 868965ea79SLois Curfman McInnes for (j=0; j<n; j++) { 875ef9f2a5SBarry Smith if (idxn[j] < 0) continue; 88273d9f13SBarry Smith if (idxn[j] >= mat->N) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large"); 8939b7565bSBarry Smith ierr = MatSetValues(A->A,1,&row,1,&idxn[j],v+i+j*m,addv);CHKERRQ(ierr); 9039b7565bSBarry Smith } 918965ea79SLois Curfman McInnes } 923a40ed3dSBarry Smith } else { 933782ba37SSatish Balay if (!A->donotstash) { 9439b7565bSBarry Smith if (roworiented) { 958798bf22SSatish Balay ierr = MatStashValuesRow_Private(&mat->stash,idxm[i],n,idxn,v+i*n);CHKERRQ(ierr); 96d36fbae8SSatish Balay } else { 978798bf22SSatish Balay ierr = MatStashValuesCol_Private(&mat->stash,idxm[i],n,idxn,v+i,m);CHKERRQ(ierr); 9839b7565bSBarry Smith } 99b49de8d1SLois Curfman McInnes } 100b49de8d1SLois Curfman McInnes } 1013782ba37SSatish Balay } 1023a40ed3dSBarry Smith PetscFunctionReturn(0); 103b49de8d1SLois Curfman McInnes } 104b49de8d1SLois Curfman McInnes 1054a2ae208SSatish Balay #undef __FUNCT__ 1064a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_MPIDense" 10713f74950SBarry Smith PetscErrorCode MatGetValues_MPIDense(Mat mat,PetscInt m,const PetscInt idxm[],PetscInt n,const PetscInt idxn[],PetscScalar v[]) 108b49de8d1SLois Curfman McInnes { 109b49de8d1SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 110dfbe8321SBarry Smith PetscErrorCode ierr; 11113f74950SBarry Smith PetscInt i,j,rstart = mdn->rstart,rend = mdn->rend,row; 112b49de8d1SLois Curfman McInnes 1133a40ed3dSBarry Smith PetscFunctionBegin; 114b49de8d1SLois Curfman McInnes for (i=0; i<m; i++) { 11529bbc08cSBarry Smith if (idxm[i] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative row"); 116273d9f13SBarry Smith if (idxm[i] >= mat->M) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large"); 117b49de8d1SLois Curfman McInnes if (idxm[i] >= rstart && idxm[i] < rend) { 118b49de8d1SLois Curfman McInnes row = idxm[i] - rstart; 119b49de8d1SLois Curfman McInnes for (j=0; j<n; j++) { 12029bbc08cSBarry Smith if (idxn[j] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative column"); 121273d9f13SBarry Smith if (idxn[j] >= mat->N) { 12229bbc08cSBarry Smith SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large"); 123a8c6a408SBarry Smith } 124b49de8d1SLois Curfman McInnes ierr = MatGetValues(mdn->A,1,&row,1,&idxn[j],v+i*n+j);CHKERRQ(ierr); 125b49de8d1SLois Curfman McInnes } 126a8c6a408SBarry Smith } else { 12729bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"Only local values currently supported"); 1288965ea79SLois Curfman McInnes } 1298965ea79SLois Curfman McInnes } 1303a40ed3dSBarry Smith PetscFunctionReturn(0); 1318965ea79SLois Curfman McInnes } 1328965ea79SLois Curfman McInnes 1334a2ae208SSatish Balay #undef __FUNCT__ 1344a2ae208SSatish Balay #define __FUNCT__ "MatGetArray_MPIDense" 135dfbe8321SBarry Smith PetscErrorCode MatGetArray_MPIDense(Mat A,PetscScalar *array[]) 136ff14e315SSatish Balay { 137ff14e315SSatish Balay Mat_MPIDense *a = (Mat_MPIDense*)A->data; 138dfbe8321SBarry Smith PetscErrorCode ierr; 139ff14e315SSatish Balay 1403a40ed3dSBarry Smith PetscFunctionBegin; 141ff14e315SSatish Balay ierr = MatGetArray(a->A,array);CHKERRQ(ierr); 1423a40ed3dSBarry Smith PetscFunctionReturn(0); 143ff14e315SSatish Balay } 144ff14e315SSatish Balay 1454a2ae208SSatish Balay #undef __FUNCT__ 1464a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrix_MPIDense" 14713f74950SBarry Smith static PetscErrorCode MatGetSubMatrix_MPIDense(Mat A,IS isrow,IS iscol,PetscInt cs,MatReuse scall,Mat *B) 148ca3fa75bSLois Curfman McInnes { 149ca3fa75bSLois Curfman McInnes Mat_MPIDense *mat = (Mat_MPIDense*)A->data,*newmatd; 150ca3fa75bSLois Curfman McInnes Mat_SeqDense *lmat = (Mat_SeqDense*)mat->A->data; 1516849ba73SBarry Smith PetscErrorCode ierr; 15213f74950SBarry Smith PetscInt i,j,*irow,*icol,rstart,rend,nrows,ncols,nlrows,nlcols; 15387828ca2SBarry Smith PetscScalar *av,*bv,*v = lmat->v; 154ca3fa75bSLois Curfman McInnes Mat newmat; 155ca3fa75bSLois Curfman McInnes 156ca3fa75bSLois Curfman McInnes PetscFunctionBegin; 157ca3fa75bSLois Curfman McInnes ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr); 158ca3fa75bSLois Curfman McInnes ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr); 159b9b97703SBarry Smith ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr); 160b9b97703SBarry Smith ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr); 161ca3fa75bSLois Curfman McInnes 162ca3fa75bSLois Curfman McInnes /* No parallel redistribution currently supported! Should really check each index set 1637eba5e9cSLois Curfman McInnes to comfirm that it is OK. ... Currently supports only submatrix same partitioning as 1647eba5e9cSLois Curfman McInnes original matrix! */ 165ca3fa75bSLois Curfman McInnes 166ca3fa75bSLois Curfman McInnes ierr = MatGetLocalSize(A,&nlrows,&nlcols);CHKERRQ(ierr); 1677eba5e9cSLois Curfman McInnes ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr); 168ca3fa75bSLois Curfman McInnes 169ca3fa75bSLois Curfman McInnes /* Check submatrix call */ 170ca3fa75bSLois Curfman McInnes if (scall == MAT_REUSE_MATRIX) { 17129bbc08cSBarry Smith /* SETERRQ(PETSC_ERR_ARG_SIZ,"Reused submatrix wrong size"); */ 1727eba5e9cSLois Curfman McInnes /* Really need to test rows and column sizes! */ 173ca3fa75bSLois Curfman McInnes newmat = *B; 174ca3fa75bSLois Curfman McInnes } else { 175ca3fa75bSLois Curfman McInnes /* Create and fill new matrix */ 176f69a0ea3SMatthew Knepley ierr = MatCreate(A->comm,&newmat);CHKERRQ(ierr); 177f69a0ea3SMatthew Knepley ierr = MatSetSizes(newmat,nrows,cs,PETSC_DECIDE,ncols);CHKERRQ(ierr); 178878740d9SKris Buschelman ierr = MatSetType(newmat,A->type_name);CHKERRQ(ierr); 179878740d9SKris Buschelman ierr = MatMPIDenseSetPreallocation(newmat,PETSC_NULL);CHKERRQ(ierr); 180ca3fa75bSLois Curfman McInnes } 181ca3fa75bSLois Curfman McInnes 182ca3fa75bSLois Curfman McInnes /* Now extract the data pointers and do the copy, column at a time */ 183ca3fa75bSLois Curfman McInnes newmatd = (Mat_MPIDense*)newmat->data; 184ca3fa75bSLois Curfman McInnes bv = ((Mat_SeqDense *)newmatd->A->data)->v; 185ca3fa75bSLois Curfman McInnes 186ca3fa75bSLois Curfman McInnes for (i=0; i<ncols; i++) { 187ca3fa75bSLois Curfman McInnes av = v + nlrows*icol[i]; 188ca3fa75bSLois Curfman McInnes for (j=0; j<nrows; j++) { 1897eba5e9cSLois Curfman McInnes *bv++ = av[irow[j] - rstart]; 190ca3fa75bSLois Curfman McInnes } 191ca3fa75bSLois Curfman McInnes } 192ca3fa75bSLois Curfman McInnes 193ca3fa75bSLois Curfman McInnes /* Assemble the matrices so that the correct flags are set */ 194ca3fa75bSLois Curfman McInnes ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 195ca3fa75bSLois Curfman McInnes ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 196ca3fa75bSLois Curfman McInnes 197ca3fa75bSLois Curfman McInnes /* Free work space */ 198ca3fa75bSLois Curfman McInnes ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr); 199ca3fa75bSLois Curfman McInnes ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr); 200ca3fa75bSLois Curfman McInnes *B = newmat; 201ca3fa75bSLois Curfman McInnes PetscFunctionReturn(0); 202ca3fa75bSLois Curfman McInnes } 203ca3fa75bSLois Curfman McInnes 2044a2ae208SSatish Balay #undef __FUNCT__ 2054a2ae208SSatish Balay #define __FUNCT__ "MatRestoreArray_MPIDense" 206dfbe8321SBarry Smith PetscErrorCode MatRestoreArray_MPIDense(Mat A,PetscScalar *array[]) 207ff14e315SSatish Balay { 2083a40ed3dSBarry Smith PetscFunctionBegin; 2093a40ed3dSBarry Smith PetscFunctionReturn(0); 210ff14e315SSatish Balay } 211ff14e315SSatish Balay 2124a2ae208SSatish Balay #undef __FUNCT__ 2134a2ae208SSatish Balay #define __FUNCT__ "MatAssemblyBegin_MPIDense" 214dfbe8321SBarry Smith PetscErrorCode MatAssemblyBegin_MPIDense(Mat mat,MatAssemblyType mode) 2158965ea79SLois Curfman McInnes { 21639ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 2178965ea79SLois Curfman McInnes MPI_Comm comm = mat->comm; 218dfbe8321SBarry Smith PetscErrorCode ierr; 21913f74950SBarry Smith PetscInt nstash,reallocs; 2208965ea79SLois Curfman McInnes InsertMode addv; 2218965ea79SLois Curfman McInnes 2223a40ed3dSBarry Smith PetscFunctionBegin; 2238965ea79SLois Curfman McInnes /* make sure all processors are either in INSERTMODE or ADDMODE */ 224ca161407SBarry Smith ierr = MPI_Allreduce(&mat->insertmode,&addv,1,MPI_INT,MPI_BOR,comm);CHKERRQ(ierr); 2257056b6fcSBarry Smith if (addv == (ADD_VALUES|INSERT_VALUES)) { 22629bbc08cSBarry Smith SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Cannot mix adds/inserts on different procs"); 2278965ea79SLois Curfman McInnes } 228e0fa3b82SLois Curfman McInnes mat->insertmode = addv; /* in case this processor had no cache */ 2298965ea79SLois Curfman McInnes 2308798bf22SSatish Balay ierr = MatStashScatterBegin_Private(&mat->stash,mdn->rowners);CHKERRQ(ierr); 2318798bf22SSatish Balay ierr = MatStashGetInfo_Private(&mat->stash,&nstash,&reallocs);CHKERRQ(ierr); 23263ba0a88SBarry Smith ierr = PetscLogInfo((mdn->A,"MatAssemblyBegin_MPIDense:Stash has %D entries, uses %D mallocs.\n",nstash,reallocs));CHKERRQ(ierr); 2333a40ed3dSBarry Smith PetscFunctionReturn(0); 2348965ea79SLois Curfman McInnes } 2358965ea79SLois Curfman McInnes 2364a2ae208SSatish Balay #undef __FUNCT__ 2374a2ae208SSatish Balay #define __FUNCT__ "MatAssemblyEnd_MPIDense" 238dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_MPIDense(Mat mat,MatAssemblyType mode) 2398965ea79SLois Curfman McInnes { 24039ddd567SLois Curfman McInnes Mat_MPIDense *mdn=(Mat_MPIDense*)mat->data; 2416849ba73SBarry Smith PetscErrorCode ierr; 24213f74950SBarry Smith PetscInt i,*row,*col,flg,j,rstart,ncols; 24313f74950SBarry Smith PetscMPIInt n; 24487828ca2SBarry Smith PetscScalar *val; 245e0fa3b82SLois Curfman McInnes InsertMode addv=mat->insertmode; 2468965ea79SLois Curfman McInnes 2473a40ed3dSBarry Smith PetscFunctionBegin; 2488965ea79SLois Curfman McInnes /* wait on receives */ 2497ef1d9bdSSatish Balay while (1) { 2508798bf22SSatish Balay ierr = MatStashScatterGetMesg_Private(&mat->stash,&n,&row,&col,&val,&flg);CHKERRQ(ierr); 2517ef1d9bdSSatish Balay if (!flg) break; 2528965ea79SLois Curfman McInnes 2537ef1d9bdSSatish Balay for (i=0; i<n;) { 2547ef1d9bdSSatish Balay /* Now identify the consecutive vals belonging to the same row */ 2557ef1d9bdSSatish Balay for (j=i,rstart=row[j]; j<n; j++) { if (row[j] != rstart) break; } 2567ef1d9bdSSatish Balay if (j < n) ncols = j-i; 2577ef1d9bdSSatish Balay else ncols = n-i; 2587ef1d9bdSSatish Balay /* Now assemble all these values with a single function call */ 2597ef1d9bdSSatish Balay ierr = MatSetValues_MPIDense(mat,1,row+i,ncols,col+i,val+i,addv);CHKERRQ(ierr); 2607ef1d9bdSSatish Balay i = j; 2618965ea79SLois Curfman McInnes } 2627ef1d9bdSSatish Balay } 2638798bf22SSatish Balay ierr = MatStashScatterEnd_Private(&mat->stash);CHKERRQ(ierr); 2648965ea79SLois Curfman McInnes 26539ddd567SLois Curfman McInnes ierr = MatAssemblyBegin(mdn->A,mode);CHKERRQ(ierr); 26639ddd567SLois Curfman McInnes ierr = MatAssemblyEnd(mdn->A,mode);CHKERRQ(ierr); 2678965ea79SLois Curfman McInnes 2686d4a8577SBarry Smith if (!mat->was_assembled && mode == MAT_FINAL_ASSEMBLY) { 26939ddd567SLois Curfman McInnes ierr = MatSetUpMultiply_MPIDense(mat);CHKERRQ(ierr); 2708965ea79SLois Curfman McInnes } 2713a40ed3dSBarry Smith PetscFunctionReturn(0); 2728965ea79SLois Curfman McInnes } 2738965ea79SLois Curfman McInnes 2744a2ae208SSatish Balay #undef __FUNCT__ 2754a2ae208SSatish Balay #define __FUNCT__ "MatZeroEntries_MPIDense" 276dfbe8321SBarry Smith PetscErrorCode MatZeroEntries_MPIDense(Mat A) 2778965ea79SLois Curfman McInnes { 278dfbe8321SBarry Smith PetscErrorCode ierr; 27939ddd567SLois Curfman McInnes Mat_MPIDense *l = (Mat_MPIDense*)A->data; 2803a40ed3dSBarry Smith 2813a40ed3dSBarry Smith PetscFunctionBegin; 2823a40ed3dSBarry Smith ierr = MatZeroEntries(l->A);CHKERRQ(ierr); 2833a40ed3dSBarry Smith PetscFunctionReturn(0); 2848965ea79SLois Curfman McInnes } 2858965ea79SLois Curfman McInnes 2868965ea79SLois Curfman McInnes /* the code does not do the diagonal entries correctly unless the 2878965ea79SLois Curfman McInnes matrix is square and the column and row owerships are identical. 2888965ea79SLois Curfman McInnes This is a BUG. The only way to fix it seems to be to access 2893501a2bdSLois Curfman McInnes mdn->A and mdn->B directly and not through the MatZeroRows() 2908965ea79SLois Curfman McInnes routine. 2918965ea79SLois Curfman McInnes */ 2924a2ae208SSatish Balay #undef __FUNCT__ 2934a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_MPIDense" 294*f4df32b1SMatthew Knepley PetscErrorCode MatZeroRows_MPIDense(Mat A,PetscInt N,const PetscInt rows[],PetscScalar diag) 2958965ea79SLois Curfman McInnes { 29639ddd567SLois Curfman McInnes Mat_MPIDense *l = (Mat_MPIDense*)A->data; 2976849ba73SBarry Smith PetscErrorCode ierr; 298*f4df32b1SMatthew Knepley PetscInt i,*owners = l->rowners; 29913f74950SBarry Smith PetscInt *nprocs,j,idx,nsends; 30013f74950SBarry Smith PetscInt nmax,*svalues,*starts,*owner,nrecvs; 30113f74950SBarry Smith PetscInt *rvalues,tag = A->tag,count,base,slen,*source; 30213f74950SBarry Smith PetscInt *lens,*lrows,*values; 30313f74950SBarry Smith PetscMPIInt n,imdex,rank = l->rank,size = l->size; 3048965ea79SLois Curfman McInnes MPI_Comm comm = A->comm; 3058965ea79SLois Curfman McInnes MPI_Request *send_waits,*recv_waits; 3068965ea79SLois Curfman McInnes MPI_Status recv_status,*send_status; 30735d8aa7fSBarry Smith PetscTruth found; 3088965ea79SLois Curfman McInnes 3093a40ed3dSBarry Smith PetscFunctionBegin; 3108965ea79SLois Curfman McInnes /* first count number of contributors to each processor */ 31113f74950SBarry Smith ierr = PetscMalloc(2*size*sizeof(PetscInt),&nprocs);CHKERRQ(ierr); 31213f74950SBarry Smith ierr = PetscMemzero(nprocs,2*size*sizeof(PetscInt));CHKERRQ(ierr); 31313f74950SBarry Smith ierr = PetscMalloc((N+1)*sizeof(PetscInt),&owner);CHKERRQ(ierr); /* see note*/ 3148965ea79SLois Curfman McInnes for (i=0; i<N; i++) { 3158965ea79SLois Curfman McInnes idx = rows[i]; 31635d8aa7fSBarry Smith found = PETSC_FALSE; 3178965ea79SLois Curfman McInnes for (j=0; j<size; j++) { 3188965ea79SLois Curfman McInnes if (idx >= owners[j] && idx < owners[j+1]) { 319c1dc657dSBarry Smith nprocs[2*j]++; nprocs[2*j+1] = 1; owner[i] = j; found = PETSC_TRUE; break; 3208965ea79SLois Curfman McInnes } 3218965ea79SLois Curfman McInnes } 32229bbc08cSBarry Smith if (!found) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Index out of range"); 3238965ea79SLois Curfman McInnes } 324c1dc657dSBarry Smith nsends = 0; for (i=0; i<size; i++) { nsends += nprocs[2*i+1];} 3258965ea79SLois Curfman McInnes 3268965ea79SLois Curfman McInnes /* inform other processors of number of messages and max length*/ 327c1dc657dSBarry Smith ierr = PetscMaxSum(comm,nprocs,&nmax,&nrecvs);CHKERRQ(ierr); 3288965ea79SLois Curfman McInnes 3298965ea79SLois Curfman McInnes /* post receives: */ 33013f74950SBarry Smith ierr = PetscMalloc((nrecvs+1)*(nmax+1)*sizeof(PetscInt),&rvalues);CHKERRQ(ierr); 331b0a32e0cSBarry Smith ierr = PetscMalloc((nrecvs+1)*sizeof(MPI_Request),&recv_waits);CHKERRQ(ierr); 3328965ea79SLois Curfman McInnes for (i=0; i<nrecvs; i++) { 33313f74950SBarry Smith ierr = MPI_Irecv(rvalues+nmax*i,nmax,MPIU_INT,MPI_ANY_SOURCE,tag,comm,recv_waits+i);CHKERRQ(ierr); 3348965ea79SLois Curfman McInnes } 3358965ea79SLois Curfman McInnes 3368965ea79SLois Curfman McInnes /* do sends: 3378965ea79SLois Curfman McInnes 1) starts[i] gives the starting index in svalues for stuff going to 3388965ea79SLois Curfman McInnes the ith processor 3398965ea79SLois Curfman McInnes */ 34013f74950SBarry Smith ierr = PetscMalloc((N+1)*sizeof(PetscInt),&svalues);CHKERRQ(ierr); 341b0a32e0cSBarry Smith ierr = PetscMalloc((nsends+1)*sizeof(MPI_Request),&send_waits);CHKERRQ(ierr); 34213f74950SBarry Smith ierr = PetscMalloc((size+1)*sizeof(PetscInt),&starts);CHKERRQ(ierr); 3438965ea79SLois Curfman McInnes starts[0] = 0; 344c1dc657dSBarry Smith for (i=1; i<size; i++) { starts[i] = starts[i-1] + nprocs[2*i-2];} 3458965ea79SLois Curfman McInnes for (i=0; i<N; i++) { 3468965ea79SLois Curfman McInnes svalues[starts[owner[i]]++] = rows[i]; 3478965ea79SLois Curfman McInnes } 3488965ea79SLois Curfman McInnes 3498965ea79SLois Curfman McInnes starts[0] = 0; 350c1dc657dSBarry Smith for (i=1; i<size+1; i++) { starts[i] = starts[i-1] + nprocs[2*i-2];} 3518965ea79SLois Curfman McInnes count = 0; 3528965ea79SLois Curfman McInnes for (i=0; i<size; i++) { 353c1dc657dSBarry Smith if (nprocs[2*i+1]) { 35413f74950SBarry Smith ierr = MPI_Isend(svalues+starts[i],nprocs[2*i],MPIU_INT,i,tag,comm,send_waits+count++);CHKERRQ(ierr); 3558965ea79SLois Curfman McInnes } 3568965ea79SLois Curfman McInnes } 357606d414cSSatish Balay ierr = PetscFree(starts);CHKERRQ(ierr); 3588965ea79SLois Curfman McInnes 3598965ea79SLois Curfman McInnes base = owners[rank]; 3608965ea79SLois Curfman McInnes 3618965ea79SLois Curfman McInnes /* wait on receives */ 36213f74950SBarry Smith ierr = PetscMalloc(2*(nrecvs+1)*sizeof(PetscInt),&lens);CHKERRQ(ierr); 3638965ea79SLois Curfman McInnes source = lens + nrecvs; 3648965ea79SLois Curfman McInnes count = nrecvs; slen = 0; 3658965ea79SLois Curfman McInnes while (count) { 366ca161407SBarry Smith ierr = MPI_Waitany(nrecvs,recv_waits,&imdex,&recv_status);CHKERRQ(ierr); 3678965ea79SLois Curfman McInnes /* unpack receives into our local space */ 36813f74950SBarry Smith ierr = MPI_Get_count(&recv_status,MPIU_INT,&n);CHKERRQ(ierr); 3698965ea79SLois Curfman McInnes source[imdex] = recv_status.MPI_SOURCE; 3708965ea79SLois Curfman McInnes lens[imdex] = n; 3718965ea79SLois Curfman McInnes slen += n; 3728965ea79SLois Curfman McInnes count--; 3738965ea79SLois Curfman McInnes } 374606d414cSSatish Balay ierr = PetscFree(recv_waits);CHKERRQ(ierr); 3758965ea79SLois Curfman McInnes 3768965ea79SLois Curfman McInnes /* move the data into the send scatter */ 37713f74950SBarry Smith ierr = PetscMalloc((slen+1)*sizeof(PetscInt),&lrows);CHKERRQ(ierr); 3788965ea79SLois Curfman McInnes count = 0; 3798965ea79SLois Curfman McInnes for (i=0; i<nrecvs; i++) { 3808965ea79SLois Curfman McInnes values = rvalues + i*nmax; 3818965ea79SLois Curfman McInnes for (j=0; j<lens[i]; j++) { 3828965ea79SLois Curfman McInnes lrows[count++] = values[j] - base; 3838965ea79SLois Curfman McInnes } 3848965ea79SLois Curfman McInnes } 385606d414cSSatish Balay ierr = PetscFree(rvalues);CHKERRQ(ierr); 386606d414cSSatish Balay ierr = PetscFree(lens);CHKERRQ(ierr); 387606d414cSSatish Balay ierr = PetscFree(owner);CHKERRQ(ierr); 388606d414cSSatish Balay ierr = PetscFree(nprocs);CHKERRQ(ierr); 3898965ea79SLois Curfman McInnes 3908965ea79SLois Curfman McInnes /* actually zap the local rows */ 391*f4df32b1SMatthew Knepley ierr = MatZeroRows(l->A,slen,lrows,diag);CHKERRQ(ierr); 392606d414cSSatish Balay ierr = PetscFree(lrows);CHKERRQ(ierr); 3938965ea79SLois Curfman McInnes 3948965ea79SLois Curfman McInnes /* wait on sends */ 3958965ea79SLois Curfman McInnes if (nsends) { 396b0a32e0cSBarry Smith ierr = PetscMalloc(nsends*sizeof(MPI_Status),&send_status);CHKERRQ(ierr); 397ca161407SBarry Smith ierr = MPI_Waitall(nsends,send_waits,send_status);CHKERRQ(ierr); 398606d414cSSatish Balay ierr = PetscFree(send_status);CHKERRQ(ierr); 3998965ea79SLois Curfman McInnes } 400606d414cSSatish Balay ierr = PetscFree(send_waits);CHKERRQ(ierr); 401606d414cSSatish Balay ierr = PetscFree(svalues);CHKERRQ(ierr); 4028965ea79SLois Curfman McInnes 4033a40ed3dSBarry Smith PetscFunctionReturn(0); 4048965ea79SLois Curfman McInnes } 4058965ea79SLois Curfman McInnes 4064a2ae208SSatish Balay #undef __FUNCT__ 4074a2ae208SSatish Balay #define __FUNCT__ "MatMult_MPIDense" 408dfbe8321SBarry Smith PetscErrorCode MatMult_MPIDense(Mat mat,Vec xx,Vec yy) 4098965ea79SLois Curfman McInnes { 41039ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 411dfbe8321SBarry Smith PetscErrorCode ierr; 412c456f294SBarry Smith 4133a40ed3dSBarry Smith PetscFunctionBegin; 41443a90d84SBarry Smith ierr = VecScatterBegin(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr); 41543a90d84SBarry Smith ierr = VecScatterEnd(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr); 41644cd7ae7SLois Curfman McInnes ierr = MatMult_SeqDense(mdn->A,mdn->lvec,yy);CHKERRQ(ierr); 4173a40ed3dSBarry Smith PetscFunctionReturn(0); 4188965ea79SLois Curfman McInnes } 4198965ea79SLois Curfman McInnes 4204a2ae208SSatish Balay #undef __FUNCT__ 4214a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_MPIDense" 422dfbe8321SBarry Smith PetscErrorCode MatMultAdd_MPIDense(Mat mat,Vec xx,Vec yy,Vec zz) 4238965ea79SLois Curfman McInnes { 42439ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 425dfbe8321SBarry Smith PetscErrorCode ierr; 426c456f294SBarry Smith 4273a40ed3dSBarry Smith PetscFunctionBegin; 42843a90d84SBarry Smith ierr = VecScatterBegin(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr); 42943a90d84SBarry Smith ierr = VecScatterEnd(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr); 43044cd7ae7SLois Curfman McInnes ierr = MatMultAdd_SeqDense(mdn->A,mdn->lvec,yy,zz);CHKERRQ(ierr); 4313a40ed3dSBarry Smith PetscFunctionReturn(0); 4328965ea79SLois Curfman McInnes } 4338965ea79SLois Curfman McInnes 4344a2ae208SSatish Balay #undef __FUNCT__ 4354a2ae208SSatish Balay #define __FUNCT__ "MatMultTranspose_MPIDense" 436dfbe8321SBarry Smith PetscErrorCode MatMultTranspose_MPIDense(Mat A,Vec xx,Vec yy) 437096963f5SLois Curfman McInnes { 438096963f5SLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense*)A->data; 439dfbe8321SBarry Smith PetscErrorCode ierr; 44087828ca2SBarry Smith PetscScalar zero = 0.0; 441096963f5SLois Curfman McInnes 4423a40ed3dSBarry Smith PetscFunctionBegin; 4432dcb1b2aSMatthew Knepley ierr = VecSet(yy,zero);CHKERRQ(ierr); 4447c922b88SBarry Smith ierr = MatMultTranspose_SeqDense(a->A,xx,a->lvec);CHKERRQ(ierr); 445537820f0SBarry Smith ierr = VecScatterBegin(a->lvec,yy,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr); 446537820f0SBarry Smith ierr = VecScatterEnd(a->lvec,yy,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr); 4473a40ed3dSBarry Smith PetscFunctionReturn(0); 448096963f5SLois Curfman McInnes } 449096963f5SLois Curfman McInnes 4504a2ae208SSatish Balay #undef __FUNCT__ 4514a2ae208SSatish Balay #define __FUNCT__ "MatMultTransposeAdd_MPIDense" 452dfbe8321SBarry Smith PetscErrorCode MatMultTransposeAdd_MPIDense(Mat A,Vec xx,Vec yy,Vec zz) 453096963f5SLois Curfman McInnes { 454096963f5SLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense*)A->data; 455dfbe8321SBarry Smith PetscErrorCode ierr; 456096963f5SLois Curfman McInnes 4573a40ed3dSBarry Smith PetscFunctionBegin; 4583501a2bdSLois Curfman McInnes ierr = VecCopy(yy,zz);CHKERRQ(ierr); 4597c922b88SBarry Smith ierr = MatMultTranspose_SeqDense(a->A,xx,a->lvec);CHKERRQ(ierr); 460537820f0SBarry Smith ierr = VecScatterBegin(a->lvec,zz,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr); 461537820f0SBarry Smith ierr = VecScatterEnd(a->lvec,zz,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr); 4623a40ed3dSBarry Smith PetscFunctionReturn(0); 463096963f5SLois Curfman McInnes } 464096963f5SLois Curfman McInnes 4654a2ae208SSatish Balay #undef __FUNCT__ 4664a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonal_MPIDense" 467dfbe8321SBarry Smith PetscErrorCode MatGetDiagonal_MPIDense(Mat A,Vec v) 4688965ea79SLois Curfman McInnes { 46939ddd567SLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense*)A->data; 470096963f5SLois Curfman McInnes Mat_SeqDense *aloc = (Mat_SeqDense*)a->A->data; 471dfbe8321SBarry Smith PetscErrorCode ierr; 47213f74950SBarry Smith PetscInt len,i,n,m = A->m,radd; 47387828ca2SBarry Smith PetscScalar *x,zero = 0.0; 474ed3cc1f0SBarry Smith 4753a40ed3dSBarry Smith PetscFunctionBegin; 4762dcb1b2aSMatthew Knepley ierr = VecSet(v,zero);CHKERRQ(ierr); 4771ebc52fbSHong Zhang ierr = VecGetArray(v,&x);CHKERRQ(ierr); 478096963f5SLois Curfman McInnes ierr = VecGetSize(v,&n);CHKERRQ(ierr); 479273d9f13SBarry Smith if (n != A->M) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec"); 480273d9f13SBarry Smith len = PetscMin(a->A->m,a->A->n); 4817ddc982cSLois Curfman McInnes radd = a->rstart*m; 48244cd7ae7SLois Curfman McInnes for (i=0; i<len; i++) { 483096963f5SLois Curfman McInnes x[i] = aloc->v[radd + i*m + i]; 484096963f5SLois Curfman McInnes } 4851ebc52fbSHong Zhang ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 4863a40ed3dSBarry Smith PetscFunctionReturn(0); 4878965ea79SLois Curfman McInnes } 4888965ea79SLois Curfman McInnes 4894a2ae208SSatish Balay #undef __FUNCT__ 4904a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_MPIDense" 491dfbe8321SBarry Smith PetscErrorCode MatDestroy_MPIDense(Mat mat) 4928965ea79SLois Curfman McInnes { 4933501a2bdSLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 494dfbe8321SBarry Smith PetscErrorCode ierr; 495ed3cc1f0SBarry Smith 4963a40ed3dSBarry Smith PetscFunctionBegin; 49794d884c6SBarry Smith 498aa482453SBarry Smith #if defined(PETSC_USE_LOG) 49977431f27SBarry Smith PetscLogObjectState((PetscObject)mat,"Rows=%D, Cols=%D",mat->M,mat->N); 5008965ea79SLois Curfman McInnes #endif 5018798bf22SSatish Balay ierr = MatStashDestroy_Private(&mat->stash);CHKERRQ(ierr); 502606d414cSSatish Balay ierr = PetscFree(mdn->rowners);CHKERRQ(ierr); 5033501a2bdSLois Curfman McInnes ierr = MatDestroy(mdn->A);CHKERRQ(ierr); 5043501a2bdSLois Curfman McInnes if (mdn->lvec) VecDestroy(mdn->lvec); 5053501a2bdSLois Curfman McInnes if (mdn->Mvctx) VecScatterDestroy(mdn->Mvctx); 506622d7880SLois Curfman McInnes if (mdn->factor) { 507606d414cSSatish Balay if (mdn->factor->temp) {ierr = PetscFree(mdn->factor->temp);CHKERRQ(ierr);} 508606d414cSSatish Balay if (mdn->factor->tag) {ierr = PetscFree(mdn->factor->tag);CHKERRQ(ierr);} 509606d414cSSatish Balay if (mdn->factor->pivots) {ierr = PetscFree(mdn->factor->pivots);CHKERRQ(ierr);} 510606d414cSSatish Balay ierr = PetscFree(mdn->factor);CHKERRQ(ierr); 511622d7880SLois Curfman McInnes } 512606d414cSSatish Balay ierr = PetscFree(mdn);CHKERRQ(ierr); 513901853e0SKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatGetDiagonalBlock_C","",PETSC_NULL);CHKERRQ(ierr); 514901853e0SKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMPIDenseSetPreallocation_C","",PETSC_NULL);CHKERRQ(ierr); 5153a40ed3dSBarry Smith PetscFunctionReturn(0); 5168965ea79SLois Curfman McInnes } 51739ddd567SLois Curfman McInnes 5184a2ae208SSatish Balay #undef __FUNCT__ 5194a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIDense_Binary" 5206849ba73SBarry Smith static PetscErrorCode MatView_MPIDense_Binary(Mat mat,PetscViewer viewer) 5218965ea79SLois Curfman McInnes { 52239ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 523dfbe8321SBarry Smith PetscErrorCode ierr; 5247056b6fcSBarry Smith 5253a40ed3dSBarry Smith PetscFunctionBegin; 52639ddd567SLois Curfman McInnes if (mdn->size == 1) { 52739ddd567SLois Curfman McInnes ierr = MatView(mdn->A,viewer);CHKERRQ(ierr); 5288965ea79SLois Curfman McInnes } 52929bbc08cSBarry Smith else SETERRQ(PETSC_ERR_SUP,"Only uniprocessor output supported"); 5303a40ed3dSBarry Smith PetscFunctionReturn(0); 5318965ea79SLois Curfman McInnes } 5328965ea79SLois Curfman McInnes 5334a2ae208SSatish Balay #undef __FUNCT__ 5344a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIDense_ASCIIorDraworSocket" 5356849ba73SBarry Smith static PetscErrorCode MatView_MPIDense_ASCIIorDraworSocket(Mat mat,PetscViewer viewer) 5368965ea79SLois Curfman McInnes { 53739ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 538dfbe8321SBarry Smith PetscErrorCode ierr; 53932dcc486SBarry Smith PetscMPIInt size = mdn->size,rank = mdn->rank; 540b0a32e0cSBarry Smith PetscViewerType vtype; 54132077d6dSBarry Smith PetscTruth iascii,isdraw; 542b0a32e0cSBarry Smith PetscViewer sviewer; 543f3ef73ceSBarry Smith PetscViewerFormat format; 5448965ea79SLois Curfman McInnes 5453a40ed3dSBarry Smith PetscFunctionBegin; 54632077d6dSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 547fb9695e5SSatish Balay ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr); 54832077d6dSBarry Smith if (iascii) { 549b0a32e0cSBarry Smith ierr = PetscViewerGetType(viewer,&vtype);CHKERRQ(ierr); 550b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 551456192e2SBarry Smith if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) { 5524e220ebcSLois Curfman McInnes MatInfo info; 553888f2ed8SSatish Balay ierr = MatGetInfo(mat,MAT_LOCAL,&info);CHKERRQ(ierr); 55477431f27SBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer," [%d] local rows %D nz %D nz alloced %D mem %D \n",rank,mat->m, 55577431f27SBarry Smith (PetscInt)info.nz_used,(PetscInt)info.nz_allocated,(PetscInt)info.memory);CHKERRQ(ierr); 556b0a32e0cSBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 5573501a2bdSLois Curfman McInnes ierr = VecScatterView(mdn->Mvctx,viewer);CHKERRQ(ierr); 5583a40ed3dSBarry Smith PetscFunctionReturn(0); 559fb9695e5SSatish Balay } else if (format == PETSC_VIEWER_ASCII_INFO) { 5603a40ed3dSBarry Smith PetscFunctionReturn(0); 5618965ea79SLois Curfman McInnes } 562f1af5d2fSBarry Smith } else if (isdraw) { 563b0a32e0cSBarry Smith PetscDraw draw; 564f1af5d2fSBarry Smith PetscTruth isnull; 565f1af5d2fSBarry Smith 566b0a32e0cSBarry Smith ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 567b0a32e0cSBarry Smith ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); 568f1af5d2fSBarry Smith if (isnull) PetscFunctionReturn(0); 569f1af5d2fSBarry Smith } 57077ed5343SBarry Smith 5718965ea79SLois Curfman McInnes if (size == 1) { 57239ddd567SLois Curfman McInnes ierr = MatView(mdn->A,viewer);CHKERRQ(ierr); 5733a40ed3dSBarry Smith } else { 5748965ea79SLois Curfman McInnes /* assemble the entire matrix onto first processor. */ 5758965ea79SLois Curfman McInnes Mat A; 57613f74950SBarry Smith PetscInt M = mat->M,N = mat->N,m,row,i,nz; 577ba8c8a56SBarry Smith PetscInt *cols; 578ba8c8a56SBarry Smith PetscScalar *vals; 5798965ea79SLois Curfman McInnes 580f69a0ea3SMatthew Knepley ierr = MatCreate(mat->comm,&A);CHKERRQ(ierr); 5818965ea79SLois Curfman McInnes if (!rank) { 582f69a0ea3SMatthew Knepley ierr = MatSetSizes(A,M,N,M,N);CHKERRQ(ierr); 5833a40ed3dSBarry Smith } else { 584f69a0ea3SMatthew Knepley ierr = MatSetSizes(A,0,0,M,N);CHKERRQ(ierr); 5858965ea79SLois Curfman McInnes } 586be5d1d56SKris Buschelman /* Since this is a temporary matrix, MATMPIDENSE instead of A->type_name here is probably acceptable. */ 587878740d9SKris Buschelman ierr = MatSetType(A,MATMPIDENSE);CHKERRQ(ierr); 588878740d9SKris Buschelman ierr = MatMPIDenseSetPreallocation(A,PETSC_NULL); 58952e6d16bSBarry Smith ierr = PetscLogObjectParent(mat,A);CHKERRQ(ierr); 5908965ea79SLois Curfman McInnes 59139ddd567SLois Curfman McInnes /* Copy the matrix ... This isn't the most efficient means, 59239ddd567SLois Curfman McInnes but it's quick for now */ 59351022da4SBarry Smith A->insertmode = INSERT_VALUES; 594273d9f13SBarry Smith row = mdn->rstart; m = mdn->A->m; 5958965ea79SLois Curfman McInnes for (i=0; i<m; i++) { 596ba8c8a56SBarry Smith ierr = MatGetRow_MPIDense(mat,row,&nz,&cols,&vals);CHKERRQ(ierr); 597ba8c8a56SBarry Smith ierr = MatSetValues_MPIDense(A,1,&row,nz,cols,vals,INSERT_VALUES);CHKERRQ(ierr); 598ba8c8a56SBarry Smith ierr = MatRestoreRow_MPIDense(mat,row,&nz,&cols,&vals);CHKERRQ(ierr); 59939ddd567SLois Curfman McInnes row++; 6008965ea79SLois Curfman McInnes } 6018965ea79SLois Curfman McInnes 6026d4a8577SBarry Smith ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 6036d4a8577SBarry Smith ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 604b0a32e0cSBarry Smith ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr); 605b9b97703SBarry Smith if (!rank) { 6066831982aSBarry Smith ierr = MatView(((Mat_MPIDense*)(A->data))->A,sviewer);CHKERRQ(ierr); 6078965ea79SLois Curfman McInnes } 608b0a32e0cSBarry Smith ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr); 609b0a32e0cSBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 6108965ea79SLois Curfman McInnes ierr = MatDestroy(A);CHKERRQ(ierr); 6118965ea79SLois Curfman McInnes } 6123a40ed3dSBarry Smith PetscFunctionReturn(0); 6138965ea79SLois Curfman McInnes } 6148965ea79SLois Curfman McInnes 6154a2ae208SSatish Balay #undef __FUNCT__ 6164a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIDense" 617dfbe8321SBarry Smith PetscErrorCode MatView_MPIDense(Mat mat,PetscViewer viewer) 6188965ea79SLois Curfman McInnes { 619dfbe8321SBarry Smith PetscErrorCode ierr; 62032077d6dSBarry Smith PetscTruth iascii,isbinary,isdraw,issocket; 6218965ea79SLois Curfman McInnes 622433994e6SBarry Smith PetscFunctionBegin; 6230f5bd95cSBarry Smith 62432077d6dSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 625fb9695e5SSatish Balay ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);CHKERRQ(ierr); 626b0a32e0cSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_SOCKET,&issocket);CHKERRQ(ierr); 627fb9695e5SSatish Balay ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr); 6280f5bd95cSBarry Smith 62932077d6dSBarry Smith if (iascii || issocket || isdraw) { 630f1af5d2fSBarry Smith ierr = MatView_MPIDense_ASCIIorDraworSocket(mat,viewer);CHKERRQ(ierr); 6310f5bd95cSBarry Smith } else if (isbinary) { 6323a40ed3dSBarry Smith ierr = MatView_MPIDense_Binary(mat,viewer);CHKERRQ(ierr); 6335cd90555SBarry Smith } else { 634958c9bccSBarry Smith SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported by MPI dense matrix",((PetscObject)viewer)->type_name); 6358965ea79SLois Curfman McInnes } 6363a40ed3dSBarry Smith PetscFunctionReturn(0); 6378965ea79SLois Curfman McInnes } 6388965ea79SLois Curfman McInnes 6394a2ae208SSatish Balay #undef __FUNCT__ 6404a2ae208SSatish Balay #define __FUNCT__ "MatGetInfo_MPIDense" 641dfbe8321SBarry Smith PetscErrorCode MatGetInfo_MPIDense(Mat A,MatInfoType flag,MatInfo *info) 6428965ea79SLois Curfman McInnes { 6433501a2bdSLois Curfman McInnes Mat_MPIDense *mat = (Mat_MPIDense*)A->data; 6443501a2bdSLois Curfman McInnes Mat mdn = mat->A; 645dfbe8321SBarry Smith PetscErrorCode ierr; 646329f5518SBarry Smith PetscReal isend[5],irecv[5]; 6478965ea79SLois Curfman McInnes 6483a40ed3dSBarry Smith PetscFunctionBegin; 649273d9f13SBarry Smith info->rows_global = (double)A->M; 650273d9f13SBarry Smith info->columns_global = (double)A->N; 651273d9f13SBarry Smith info->rows_local = (double)A->m; 652273d9f13SBarry Smith info->columns_local = (double)A->N; 6534e220ebcSLois Curfman McInnes info->block_size = 1.0; 6544e220ebcSLois Curfman McInnes ierr = MatGetInfo(mdn,MAT_LOCAL,info);CHKERRQ(ierr); 6554e220ebcSLois Curfman McInnes isend[0] = info->nz_used; isend[1] = info->nz_allocated; isend[2] = info->nz_unneeded; 6564e220ebcSLois Curfman McInnes isend[3] = info->memory; isend[4] = info->mallocs; 6578965ea79SLois Curfman McInnes if (flag == MAT_LOCAL) { 6584e220ebcSLois Curfman McInnes info->nz_used = isend[0]; 6594e220ebcSLois Curfman McInnes info->nz_allocated = isend[1]; 6604e220ebcSLois Curfman McInnes info->nz_unneeded = isend[2]; 6614e220ebcSLois Curfman McInnes info->memory = isend[3]; 6624e220ebcSLois Curfman McInnes info->mallocs = isend[4]; 6638965ea79SLois Curfman McInnes } else if (flag == MAT_GLOBAL_MAX) { 664d7d1e502SBarry Smith ierr = MPI_Allreduce(isend,irecv,5,MPIU_REAL,MPI_MAX,A->comm);CHKERRQ(ierr); 6654e220ebcSLois Curfman McInnes info->nz_used = irecv[0]; 6664e220ebcSLois Curfman McInnes info->nz_allocated = irecv[1]; 6674e220ebcSLois Curfman McInnes info->nz_unneeded = irecv[2]; 6684e220ebcSLois Curfman McInnes info->memory = irecv[3]; 6694e220ebcSLois Curfman McInnes info->mallocs = irecv[4]; 6708965ea79SLois Curfman McInnes } else if (flag == MAT_GLOBAL_SUM) { 671d7d1e502SBarry Smith ierr = MPI_Allreduce(isend,irecv,5,MPIU_REAL,MPI_SUM,A->comm);CHKERRQ(ierr); 6724e220ebcSLois Curfman McInnes info->nz_used = irecv[0]; 6734e220ebcSLois Curfman McInnes info->nz_allocated = irecv[1]; 6744e220ebcSLois Curfman McInnes info->nz_unneeded = irecv[2]; 6754e220ebcSLois Curfman McInnes info->memory = irecv[3]; 6764e220ebcSLois Curfman McInnes info->mallocs = irecv[4]; 6778965ea79SLois Curfman McInnes } 6784e220ebcSLois Curfman McInnes info->fill_ratio_given = 0; /* no parallel LU/ILU/Cholesky */ 6794e220ebcSLois Curfman McInnes info->fill_ratio_needed = 0; 6804e220ebcSLois Curfman McInnes info->factor_mallocs = 0; 6813a40ed3dSBarry Smith PetscFunctionReturn(0); 6828965ea79SLois Curfman McInnes } 6838965ea79SLois Curfman McInnes 6844a2ae208SSatish Balay #undef __FUNCT__ 6854a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_MPIDense" 686dfbe8321SBarry Smith PetscErrorCode MatSetOption_MPIDense(Mat A,MatOption op) 6878965ea79SLois Curfman McInnes { 68839ddd567SLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense*)A->data; 689dfbe8321SBarry Smith PetscErrorCode ierr; 6908965ea79SLois Curfman McInnes 6913a40ed3dSBarry Smith PetscFunctionBegin; 69212c028f9SKris Buschelman switch (op) { 69312c028f9SKris Buschelman case MAT_NO_NEW_NONZERO_LOCATIONS: 69412c028f9SKris Buschelman case MAT_YES_NEW_NONZERO_LOCATIONS: 69512c028f9SKris Buschelman case MAT_NEW_NONZERO_LOCATION_ERR: 69612c028f9SKris Buschelman case MAT_NEW_NONZERO_ALLOCATION_ERR: 69712c028f9SKris Buschelman case MAT_COLUMNS_SORTED: 69812c028f9SKris Buschelman case MAT_COLUMNS_UNSORTED: 699273d9f13SBarry Smith ierr = MatSetOption(a->A,op);CHKERRQ(ierr); 70012c028f9SKris Buschelman break; 70112c028f9SKris Buschelman case MAT_ROW_ORIENTED: 702273d9f13SBarry Smith a->roworiented = PETSC_TRUE; 703273d9f13SBarry Smith ierr = MatSetOption(a->A,op);CHKERRQ(ierr); 70412c028f9SKris Buschelman break; 70512c028f9SKris Buschelman case MAT_ROWS_SORTED: 70612c028f9SKris Buschelman case MAT_ROWS_UNSORTED: 70712c028f9SKris Buschelman case MAT_YES_NEW_DIAGONALS: 70812c028f9SKris Buschelman case MAT_USE_HASH_TABLE: 70963ba0a88SBarry Smith ierr = PetscLogInfo((A,"MatSetOption_MPIDense:Option ignored\n"));CHKERRQ(ierr); 71012c028f9SKris Buschelman break; 71112c028f9SKris Buschelman case MAT_COLUMN_ORIENTED: 712273d9f13SBarry Smith a->roworiented = PETSC_FALSE; 713273d9f13SBarry Smith ierr = MatSetOption(a->A,op);CHKERRQ(ierr); 71412c028f9SKris Buschelman break; 71512c028f9SKris Buschelman case MAT_IGNORE_OFF_PROC_ENTRIES: 716273d9f13SBarry Smith a->donotstash = PETSC_TRUE; 71712c028f9SKris Buschelman break; 71812c028f9SKris Buschelman case MAT_NO_NEW_DIAGONALS: 71929bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"MAT_NO_NEW_DIAGONALS"); 72077e54ba9SKris Buschelman case MAT_SYMMETRIC: 72177e54ba9SKris Buschelman case MAT_STRUCTURALLY_SYMMETRIC: 7229a4540c5SBarry Smith case MAT_NOT_SYMMETRIC: 7239a4540c5SBarry Smith case MAT_NOT_STRUCTURALLY_SYMMETRIC: 7249a4540c5SBarry Smith case MAT_HERMITIAN: 7259a4540c5SBarry Smith case MAT_NOT_HERMITIAN: 7269a4540c5SBarry Smith case MAT_SYMMETRY_ETERNAL: 7279a4540c5SBarry Smith case MAT_NOT_SYMMETRY_ETERNAL: 72877e54ba9SKris Buschelman break; 72912c028f9SKris Buschelman default: 73029bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"unknown option"); 7313a40ed3dSBarry Smith } 7323a40ed3dSBarry Smith PetscFunctionReturn(0); 7338965ea79SLois Curfman McInnes } 7348965ea79SLois Curfman McInnes 7358965ea79SLois Curfman McInnes 7364a2ae208SSatish Balay #undef __FUNCT__ 7374a2ae208SSatish Balay #define __FUNCT__ "MatDiagonalScale_MPIDense" 738dfbe8321SBarry Smith PetscErrorCode MatDiagonalScale_MPIDense(Mat A,Vec ll,Vec rr) 7395b2fa520SLois Curfman McInnes { 7405b2fa520SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)A->data; 7415b2fa520SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)mdn->A->data; 74287828ca2SBarry Smith PetscScalar *l,*r,x,*v; 743dfbe8321SBarry Smith PetscErrorCode ierr; 74413f74950SBarry Smith PetscInt i,j,s2a,s3a,s2,s3,m=mdn->A->m,n=mdn->A->n; 7455b2fa520SLois Curfman McInnes 7465b2fa520SLois Curfman McInnes PetscFunctionBegin; 74772d926a5SLois Curfman McInnes ierr = MatGetLocalSize(A,&s2,&s3);CHKERRQ(ierr); 7485b2fa520SLois Curfman McInnes if (ll) { 74972d926a5SLois Curfman McInnes ierr = VecGetLocalSize(ll,&s2a);CHKERRQ(ierr); 75029bbc08cSBarry Smith if (s2a != s2) SETERRQ(PETSC_ERR_ARG_SIZ,"Left scaling vector non-conforming local size"); 7511ebc52fbSHong Zhang ierr = VecGetArray(ll,&l);CHKERRQ(ierr); 7525b2fa520SLois Curfman McInnes for (i=0; i<m; i++) { 7535b2fa520SLois Curfman McInnes x = l[i]; 7545b2fa520SLois Curfman McInnes v = mat->v + i; 7555b2fa520SLois Curfman McInnes for (j=0; j<n; j++) { (*v) *= x; v+= m;} 7565b2fa520SLois Curfman McInnes } 7571ebc52fbSHong Zhang ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr); 758efee365bSSatish Balay ierr = PetscLogFlops(n*m);CHKERRQ(ierr); 7595b2fa520SLois Curfman McInnes } 7605b2fa520SLois Curfman McInnes if (rr) { 76172d926a5SLois Curfman McInnes ierr = VecGetSize(rr,&s3a);CHKERRQ(ierr); 76229bbc08cSBarry Smith if (s3a != s3) SETERRQ(PETSC_ERR_ARG_SIZ,"Right scaling vec non-conforming local size"); 7635b2fa520SLois Curfman McInnes ierr = VecScatterBegin(rr,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr); 7645b2fa520SLois Curfman McInnes ierr = VecScatterEnd(rr,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr); 7651ebc52fbSHong Zhang ierr = VecGetArray(mdn->lvec,&r);CHKERRQ(ierr); 7665b2fa520SLois Curfman McInnes for (i=0; i<n; i++) { 7675b2fa520SLois Curfman McInnes x = r[i]; 7685b2fa520SLois Curfman McInnes v = mat->v + i*m; 7695b2fa520SLois Curfman McInnes for (j=0; j<m; j++) { (*v++) *= x;} 7705b2fa520SLois Curfman McInnes } 7711ebc52fbSHong Zhang ierr = VecRestoreArray(mdn->lvec,&r);CHKERRQ(ierr); 772efee365bSSatish Balay ierr = PetscLogFlops(n*m);CHKERRQ(ierr); 7735b2fa520SLois Curfman McInnes } 7745b2fa520SLois Curfman McInnes PetscFunctionReturn(0); 7755b2fa520SLois Curfman McInnes } 7765b2fa520SLois Curfman McInnes 7774a2ae208SSatish Balay #undef __FUNCT__ 7784a2ae208SSatish Balay #define __FUNCT__ "MatNorm_MPIDense" 779dfbe8321SBarry Smith PetscErrorCode MatNorm_MPIDense(Mat A,NormType type,PetscReal *nrm) 780096963f5SLois Curfman McInnes { 7813501a2bdSLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)A->data; 7823501a2bdSLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)mdn->A->data; 783dfbe8321SBarry Smith PetscErrorCode ierr; 78413f74950SBarry Smith PetscInt i,j; 785329f5518SBarry Smith PetscReal sum = 0.0; 78687828ca2SBarry Smith PetscScalar *v = mat->v; 7873501a2bdSLois Curfman McInnes 7883a40ed3dSBarry Smith PetscFunctionBegin; 7893501a2bdSLois Curfman McInnes if (mdn->size == 1) { 790064f8208SBarry Smith ierr = MatNorm(mdn->A,type,nrm);CHKERRQ(ierr); 7913501a2bdSLois Curfman McInnes } else { 7923501a2bdSLois Curfman McInnes if (type == NORM_FROBENIUS) { 793273d9f13SBarry Smith for (i=0; i<mdn->A->n*mdn->A->m; i++) { 794aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 795329f5518SBarry Smith sum += PetscRealPart(PetscConj(*v)*(*v)); v++; 7963501a2bdSLois Curfman McInnes #else 7973501a2bdSLois Curfman McInnes sum += (*v)*(*v); v++; 7983501a2bdSLois Curfman McInnes #endif 7993501a2bdSLois Curfman McInnes } 800064f8208SBarry Smith ierr = MPI_Allreduce(&sum,nrm,1,MPIU_REAL,MPI_SUM,A->comm);CHKERRQ(ierr); 801064f8208SBarry Smith *nrm = sqrt(*nrm); 802efee365bSSatish Balay ierr = PetscLogFlops(2*mdn->A->n*mdn->A->m);CHKERRQ(ierr); 8033a40ed3dSBarry Smith } else if (type == NORM_1) { 804329f5518SBarry Smith PetscReal *tmp,*tmp2; 805b0a32e0cSBarry Smith ierr = PetscMalloc(2*A->N*sizeof(PetscReal),&tmp);CHKERRQ(ierr); 806273d9f13SBarry Smith tmp2 = tmp + A->N; 807273d9f13SBarry Smith ierr = PetscMemzero(tmp,2*A->N*sizeof(PetscReal));CHKERRQ(ierr); 808064f8208SBarry Smith *nrm = 0.0; 8093501a2bdSLois Curfman McInnes v = mat->v; 810273d9f13SBarry Smith for (j=0; j<mdn->A->n; j++) { 811273d9f13SBarry Smith for (i=0; i<mdn->A->m; i++) { 81267e560aaSBarry Smith tmp[j] += PetscAbsScalar(*v); v++; 8133501a2bdSLois Curfman McInnes } 8143501a2bdSLois Curfman McInnes } 815d7d1e502SBarry Smith ierr = MPI_Allreduce(tmp,tmp2,A->N,MPIU_REAL,MPI_SUM,A->comm);CHKERRQ(ierr); 816273d9f13SBarry Smith for (j=0; j<A->N; j++) { 817064f8208SBarry Smith if (tmp2[j] > *nrm) *nrm = tmp2[j]; 8183501a2bdSLois Curfman McInnes } 819606d414cSSatish Balay ierr = PetscFree(tmp);CHKERRQ(ierr); 820efee365bSSatish Balay ierr = PetscLogFlops(A->n*A->m);CHKERRQ(ierr); 8213a40ed3dSBarry Smith } else if (type == NORM_INFINITY) { /* max row norm */ 822329f5518SBarry Smith PetscReal ntemp; 8233501a2bdSLois Curfman McInnes ierr = MatNorm(mdn->A,type,&ntemp);CHKERRQ(ierr); 824064f8208SBarry Smith ierr = MPI_Allreduce(&ntemp,nrm,1,MPIU_REAL,MPI_MAX,A->comm);CHKERRQ(ierr); 8253a40ed3dSBarry Smith } else { 82629bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"No support for two norm"); 8273501a2bdSLois Curfman McInnes } 8283501a2bdSLois Curfman McInnes } 8293a40ed3dSBarry Smith PetscFunctionReturn(0); 8303501a2bdSLois Curfman McInnes } 8313501a2bdSLois Curfman McInnes 8324a2ae208SSatish Balay #undef __FUNCT__ 8334a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_MPIDense" 834dfbe8321SBarry Smith PetscErrorCode MatTranspose_MPIDense(Mat A,Mat *matout) 8353501a2bdSLois Curfman McInnes { 8363501a2bdSLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense*)A->data; 8373501a2bdSLois Curfman McInnes Mat_SeqDense *Aloc = (Mat_SeqDense*)a->A->data; 8383501a2bdSLois Curfman McInnes Mat B; 83913f74950SBarry Smith PetscInt M = A->M,N = A->N,m,n,*rwork,rstart = a->rstart; 8406849ba73SBarry Smith PetscErrorCode ierr; 84113f74950SBarry Smith PetscInt j,i; 84287828ca2SBarry Smith PetscScalar *v; 8433501a2bdSLois Curfman McInnes 8443a40ed3dSBarry Smith PetscFunctionBegin; 8457c922b88SBarry Smith if (!matout && M != N) { 84629bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"Supports square matrix only in-place"); 8477056b6fcSBarry Smith } 848f69a0ea3SMatthew Knepley ierr = MatCreate(A->comm,&B);CHKERRQ(ierr); 849f69a0ea3SMatthew Knepley ierr = MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,N,M);CHKERRQ(ierr); 850878740d9SKris Buschelman ierr = MatSetType(B,A->type_name);CHKERRQ(ierr); 851878740d9SKris Buschelman ierr = MatMPIDenseSetPreallocation(B,PETSC_NULL);CHKERRQ(ierr); 8523501a2bdSLois Curfman McInnes 853273d9f13SBarry Smith m = a->A->m; n = a->A->n; v = Aloc->v; 85413f74950SBarry Smith ierr = PetscMalloc(n*sizeof(PetscInt),&rwork);CHKERRQ(ierr); 8553501a2bdSLois Curfman McInnes for (j=0; j<n; j++) { 8563501a2bdSLois Curfman McInnes for (i=0; i<m; i++) rwork[i] = rstart + i; 8573501a2bdSLois Curfman McInnes ierr = MatSetValues(B,1,&j,m,rwork,v,INSERT_VALUES);CHKERRQ(ierr); 8583501a2bdSLois Curfman McInnes v += m; 8593501a2bdSLois Curfman McInnes } 860606d414cSSatish Balay ierr = PetscFree(rwork);CHKERRQ(ierr); 8616d4a8577SBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 8626d4a8577SBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 8637c922b88SBarry Smith if (matout) { 8643501a2bdSLois Curfman McInnes *matout = B; 8653501a2bdSLois Curfman McInnes } else { 866273d9f13SBarry Smith ierr = MatHeaderCopy(A,B);CHKERRQ(ierr); 8673501a2bdSLois Curfman McInnes } 8683a40ed3dSBarry Smith PetscFunctionReturn(0); 869096963f5SLois Curfman McInnes } 870096963f5SLois Curfman McInnes 871d9eff348SSatish Balay #include "petscblaslapack.h" 8724a2ae208SSatish Balay #undef __FUNCT__ 8734a2ae208SSatish Balay #define __FUNCT__ "MatScale_MPIDense" 874*f4df32b1SMatthew Knepley PetscErrorCode MatScale_MPIDense(Mat inA,PetscScalar alpha) 87544cd7ae7SLois Curfman McInnes { 87644cd7ae7SLois Curfman McInnes Mat_MPIDense *A = (Mat_MPIDense*)inA->data; 87744cd7ae7SLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense*)A->A->data; 878*f4df32b1SMatthew Knepley PetscScalar oalpha = alpha; 8794ce68768SBarry Smith PetscBLASInt one = 1,nz = (PetscBLASInt)inA->m*inA->N; 880efee365bSSatish Balay PetscErrorCode ierr; 88144cd7ae7SLois Curfman McInnes 8823a40ed3dSBarry Smith PetscFunctionBegin; 883*f4df32b1SMatthew Knepley BLASscal_(&nz,&oalpha,a->v,&one); 884efee365bSSatish Balay ierr = PetscLogFlops(nz);CHKERRQ(ierr); 8853a40ed3dSBarry Smith PetscFunctionReturn(0); 88644cd7ae7SLois Curfman McInnes } 88744cd7ae7SLois Curfman McInnes 8886849ba73SBarry Smith static PetscErrorCode MatDuplicate_MPIDense(Mat,MatDuplicateOption,Mat *); 8898965ea79SLois Curfman McInnes 8904a2ae208SSatish Balay #undef __FUNCT__ 8914a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation_MPIDense" 892dfbe8321SBarry Smith PetscErrorCode MatSetUpPreallocation_MPIDense(Mat A) 893273d9f13SBarry Smith { 894dfbe8321SBarry Smith PetscErrorCode ierr; 895273d9f13SBarry Smith 896273d9f13SBarry Smith PetscFunctionBegin; 897273d9f13SBarry Smith ierr = MatMPIDenseSetPreallocation(A,0);CHKERRQ(ierr); 898273d9f13SBarry Smith PetscFunctionReturn(0); 899273d9f13SBarry Smith } 900273d9f13SBarry Smith 9018965ea79SLois Curfman McInnes /* -------------------------------------------------------------------*/ 90209dc0095SBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_MPIDense, 90309dc0095SBarry Smith MatGetRow_MPIDense, 90409dc0095SBarry Smith MatRestoreRow_MPIDense, 90509dc0095SBarry Smith MatMult_MPIDense, 90697304618SKris Buschelman /* 4*/ MatMultAdd_MPIDense, 9077c922b88SBarry Smith MatMultTranspose_MPIDense, 9087c922b88SBarry Smith MatMultTransposeAdd_MPIDense, 9098965ea79SLois Curfman McInnes 0, 91009dc0095SBarry Smith 0, 91109dc0095SBarry Smith 0, 91297304618SKris Buschelman /*10*/ 0, 91309dc0095SBarry Smith 0, 91409dc0095SBarry Smith 0, 91509dc0095SBarry Smith 0, 91609dc0095SBarry Smith MatTranspose_MPIDense, 91797304618SKris Buschelman /*15*/ MatGetInfo_MPIDense, 91897304618SKris Buschelman 0, 91909dc0095SBarry Smith MatGetDiagonal_MPIDense, 9205b2fa520SLois Curfman McInnes MatDiagonalScale_MPIDense, 92109dc0095SBarry Smith MatNorm_MPIDense, 92297304618SKris Buschelman /*20*/ MatAssemblyBegin_MPIDense, 92309dc0095SBarry Smith MatAssemblyEnd_MPIDense, 92409dc0095SBarry Smith 0, 92509dc0095SBarry Smith MatSetOption_MPIDense, 92609dc0095SBarry Smith MatZeroEntries_MPIDense, 92797304618SKris Buschelman /*25*/ MatZeroRows_MPIDense, 92809dc0095SBarry Smith 0, 92909dc0095SBarry Smith 0, 93009dc0095SBarry Smith 0, 93109dc0095SBarry Smith 0, 93297304618SKris Buschelman /*30*/ MatSetUpPreallocation_MPIDense, 933273d9f13SBarry Smith 0, 93409dc0095SBarry Smith 0, 93509dc0095SBarry Smith MatGetArray_MPIDense, 93609dc0095SBarry Smith MatRestoreArray_MPIDense, 93797304618SKris Buschelman /*35*/ MatDuplicate_MPIDense, 93809dc0095SBarry Smith 0, 93909dc0095SBarry Smith 0, 94009dc0095SBarry Smith 0, 94109dc0095SBarry Smith 0, 94297304618SKris Buschelman /*40*/ 0, 9432ce60cd0SSatish Balay MatGetSubMatrices_MPIDense, 94409dc0095SBarry Smith 0, 94509dc0095SBarry Smith MatGetValues_MPIDense, 94609dc0095SBarry Smith 0, 94797304618SKris Buschelman /*45*/ 0, 94809dc0095SBarry Smith MatScale_MPIDense, 94909dc0095SBarry Smith 0, 95009dc0095SBarry Smith 0, 95109dc0095SBarry Smith 0, 952521d7252SBarry Smith /*50*/ 0, 95309dc0095SBarry Smith 0, 95409dc0095SBarry Smith 0, 95509dc0095SBarry Smith 0, 95609dc0095SBarry Smith 0, 95797304618SKris Buschelman /*55*/ 0, 95809dc0095SBarry Smith 0, 95909dc0095SBarry Smith 0, 96009dc0095SBarry Smith 0, 96109dc0095SBarry Smith 0, 96297304618SKris Buschelman /*60*/ MatGetSubMatrix_MPIDense, 963b9b97703SBarry Smith MatDestroy_MPIDense, 964b9b97703SBarry Smith MatView_MPIDense, 96597304618SKris Buschelman MatGetPetscMaps_Petsc, 96697304618SKris Buschelman 0, 96797304618SKris Buschelman /*65*/ 0, 96897304618SKris Buschelman 0, 96997304618SKris Buschelman 0, 97097304618SKris Buschelman 0, 97197304618SKris Buschelman 0, 97297304618SKris Buschelman /*70*/ 0, 97397304618SKris Buschelman 0, 97497304618SKris Buschelman 0, 97597304618SKris Buschelman 0, 97697304618SKris Buschelman 0, 97797304618SKris Buschelman /*75*/ 0, 97897304618SKris Buschelman 0, 97997304618SKris Buschelman 0, 98097304618SKris Buschelman 0, 98197304618SKris Buschelman 0, 98297304618SKris Buschelman /*80*/ 0, 98397304618SKris Buschelman 0, 98497304618SKris Buschelman 0, 98597304618SKris Buschelman 0, 986865e5f61SKris Buschelman /*84*/ MatLoad_MPIDense, 987865e5f61SKris Buschelman 0, 988865e5f61SKris Buschelman 0, 989865e5f61SKris Buschelman 0, 990865e5f61SKris Buschelman 0, 991865e5f61SKris Buschelman 0, 992865e5f61SKris Buschelman /*90*/ 0, 993865e5f61SKris Buschelman 0, 994865e5f61SKris Buschelman 0, 995865e5f61SKris Buschelman 0, 996865e5f61SKris Buschelman 0, 997865e5f61SKris Buschelman /*95*/ 0, 998865e5f61SKris Buschelman 0, 999865e5f61SKris Buschelman 0, 1000865e5f61SKris Buschelman 0}; 10018965ea79SLois Curfman McInnes 1002273d9f13SBarry Smith EXTERN_C_BEGIN 10034a2ae208SSatish Balay #undef __FUNCT__ 1004a23d5eceSKris Buschelman #define __FUNCT__ "MatMPIDenseSetPreallocation_MPIDense" 1005be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatMPIDenseSetPreallocation_MPIDense(Mat mat,PetscScalar *data) 1006a23d5eceSKris Buschelman { 1007a23d5eceSKris Buschelman Mat_MPIDense *a; 1008dfbe8321SBarry Smith PetscErrorCode ierr; 1009a23d5eceSKris Buschelman 1010a23d5eceSKris Buschelman PetscFunctionBegin; 1011a23d5eceSKris Buschelman mat->preallocated = PETSC_TRUE; 1012a23d5eceSKris Buschelman /* Note: For now, when data is specified above, this assumes the user correctly 1013a23d5eceSKris Buschelman allocates the local dense storage space. We should add error checking. */ 1014a23d5eceSKris Buschelman 1015a23d5eceSKris Buschelman a = (Mat_MPIDense*)mat->data; 1016f69a0ea3SMatthew Knepley ierr = MatCreate(PETSC_COMM_SELF,&a->A);CHKERRQ(ierr); 1017f69a0ea3SMatthew Knepley ierr = MatSetSizes(a->A,mat->m,mat->N,mat->m,mat->N);CHKERRQ(ierr); 10185c5985e7SKris Buschelman ierr = MatSetType(a->A,MATSEQDENSE);CHKERRQ(ierr); 10195c5985e7SKris Buschelman ierr = MatSeqDenseSetPreallocation(a->A,data);CHKERRQ(ierr); 102052e6d16bSBarry Smith ierr = PetscLogObjectParent(mat,a->A);CHKERRQ(ierr); 1021a23d5eceSKris Buschelman PetscFunctionReturn(0); 1022a23d5eceSKris Buschelman } 1023a23d5eceSKris Buschelman EXTERN_C_END 1024a23d5eceSKris Buschelman 10250bad9183SKris Buschelman /*MC 1026fafad747SKris Buschelman MATMPIDENSE - MATMPIDENSE = "mpidense" - A matrix type to be used for distributed dense matrices. 10270bad9183SKris Buschelman 10280bad9183SKris Buschelman Options Database Keys: 10290bad9183SKris Buschelman . -mat_type mpidense - sets the matrix type to "mpidense" during a call to MatSetFromOptions() 10300bad9183SKris Buschelman 10310bad9183SKris Buschelman Level: beginner 10320bad9183SKris Buschelman 10330bad9183SKris Buschelman .seealso: MatCreateMPIDense 10340bad9183SKris Buschelman M*/ 10350bad9183SKris Buschelman 1036a23d5eceSKris Buschelman EXTERN_C_BEGIN 1037a23d5eceSKris Buschelman #undef __FUNCT__ 10384a2ae208SSatish Balay #define __FUNCT__ "MatCreate_MPIDense" 1039be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_MPIDense(Mat mat) 1040273d9f13SBarry Smith { 1041273d9f13SBarry Smith Mat_MPIDense *a; 1042dfbe8321SBarry Smith PetscErrorCode ierr; 104313f74950SBarry Smith PetscInt i; 1044273d9f13SBarry Smith 1045273d9f13SBarry Smith PetscFunctionBegin; 1046b0a32e0cSBarry Smith ierr = PetscNew(Mat_MPIDense,&a);CHKERRQ(ierr); 1047b0a32e0cSBarry Smith mat->data = (void*)a; 1048273d9f13SBarry Smith ierr = PetscMemcpy(mat->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 1049273d9f13SBarry Smith mat->factor = 0; 1050273d9f13SBarry Smith mat->mapping = 0; 1051273d9f13SBarry Smith 1052273d9f13SBarry Smith a->factor = 0; 1053273d9f13SBarry Smith mat->insertmode = NOT_SET_VALUES; 1054273d9f13SBarry Smith ierr = MPI_Comm_rank(mat->comm,&a->rank);CHKERRQ(ierr); 1055273d9f13SBarry Smith ierr = MPI_Comm_size(mat->comm,&a->size);CHKERRQ(ierr); 1056273d9f13SBarry Smith 1057273d9f13SBarry Smith ierr = PetscSplitOwnership(mat->comm,&mat->m,&mat->M);CHKERRQ(ierr); 1058273d9f13SBarry Smith ierr = PetscSplitOwnership(mat->comm,&mat->n,&mat->N);CHKERRQ(ierr); 1059273d9f13SBarry Smith a->nvec = mat->n; 1060273d9f13SBarry Smith 1061273d9f13SBarry Smith /* the information in the maps duplicates the information computed below, eventually 1062273d9f13SBarry Smith we should remove the duplicate information that is not contained in the maps */ 10638a124369SBarry Smith ierr = PetscMapCreateMPI(mat->comm,mat->m,mat->M,&mat->rmap);CHKERRQ(ierr); 10648a124369SBarry Smith ierr = PetscMapCreateMPI(mat->comm,mat->n,mat->N,&mat->cmap);CHKERRQ(ierr); 1065273d9f13SBarry Smith 1066273d9f13SBarry Smith /* build local table of row and column ownerships */ 106713f74950SBarry Smith ierr = PetscMalloc(2*(a->size+2)*sizeof(PetscInt),&a->rowners);CHKERRQ(ierr); 1068273d9f13SBarry Smith a->cowners = a->rowners + a->size + 1; 106952e6d16bSBarry Smith ierr = PetscLogObjectMemory(mat,2*(a->size+2)*sizeof(PetscInt)+sizeof(struct _p_Mat)+sizeof(Mat_MPIDense));CHKERRQ(ierr); 107013f74950SBarry Smith ierr = MPI_Allgather(&mat->m,1,MPIU_INT,a->rowners+1,1,MPIU_INT,mat->comm);CHKERRQ(ierr); 1071273d9f13SBarry Smith a->rowners[0] = 0; 1072273d9f13SBarry Smith for (i=2; i<=a->size; i++) { 1073273d9f13SBarry Smith a->rowners[i] += a->rowners[i-1]; 1074273d9f13SBarry Smith } 1075273d9f13SBarry Smith a->rstart = a->rowners[a->rank]; 1076273d9f13SBarry Smith a->rend = a->rowners[a->rank+1]; 107713f74950SBarry Smith ierr = MPI_Allgather(&mat->n,1,MPIU_INT,a->cowners+1,1,MPIU_INT,mat->comm);CHKERRQ(ierr); 1078273d9f13SBarry Smith a->cowners[0] = 0; 1079273d9f13SBarry Smith for (i=2; i<=a->size; i++) { 1080273d9f13SBarry Smith a->cowners[i] += a->cowners[i-1]; 1081273d9f13SBarry Smith } 1082273d9f13SBarry Smith 1083273d9f13SBarry Smith /* build cache for off array entries formed */ 1084273d9f13SBarry Smith a->donotstash = PETSC_FALSE; 1085273d9f13SBarry Smith ierr = MatStashCreate_Private(mat->comm,1,&mat->stash);CHKERRQ(ierr); 1086273d9f13SBarry Smith 1087273d9f13SBarry Smith /* stuff used for matrix vector multiply */ 1088273d9f13SBarry Smith a->lvec = 0; 1089273d9f13SBarry Smith a->Mvctx = 0; 1090273d9f13SBarry Smith a->roworiented = PETSC_TRUE; 1091273d9f13SBarry Smith 1092273d9f13SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatGetDiagonalBlock_C", 1093273d9f13SBarry Smith "MatGetDiagonalBlock_MPIDense", 1094273d9f13SBarry Smith MatGetDiagonalBlock_MPIDense);CHKERRQ(ierr); 1095a23d5eceSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMPIDenseSetPreallocation_C", 1096a23d5eceSKris Buschelman "MatMPIDenseSetPreallocation_MPIDense", 1097a23d5eceSKris Buschelman MatMPIDenseSetPreallocation_MPIDense);CHKERRQ(ierr); 1098273d9f13SBarry Smith PetscFunctionReturn(0); 1099273d9f13SBarry Smith } 1100273d9f13SBarry Smith EXTERN_C_END 1101273d9f13SBarry Smith 1102209238afSKris Buschelman /*MC 1103002d173eSKris Buschelman MATDENSE - MATDENSE = "dense" - A matrix type to be used for dense matrices. 1104209238afSKris Buschelman 1105209238afSKris Buschelman This matrix type is identical to MATSEQDENSE when constructed with a single process communicator, 1106209238afSKris Buschelman and MATMPIDENSE otherwise. 1107209238afSKris Buschelman 1108209238afSKris Buschelman Options Database Keys: 1109209238afSKris Buschelman . -mat_type dense - sets the matrix type to "dense" during a call to MatSetFromOptions() 1110209238afSKris Buschelman 1111209238afSKris Buschelman Level: beginner 1112209238afSKris Buschelman 1113209238afSKris Buschelman .seealso: MatCreateMPIDense,MATSEQDENSE,MATMPIDENSE 1114209238afSKris Buschelman M*/ 1115209238afSKris Buschelman 1116209238afSKris Buschelman EXTERN_C_BEGIN 1117209238afSKris Buschelman #undef __FUNCT__ 1118209238afSKris Buschelman #define __FUNCT__ "MatCreate_Dense" 1119be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_Dense(Mat A) 1120dfbe8321SBarry Smith { 11216849ba73SBarry Smith PetscErrorCode ierr; 112213f74950SBarry Smith PetscMPIInt size; 1123209238afSKris Buschelman 1124209238afSKris Buschelman PetscFunctionBegin; 1125209238afSKris Buschelman ierr = PetscObjectChangeTypeName((PetscObject)A,MATDENSE);CHKERRQ(ierr); 1126209238afSKris Buschelman ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr); 1127209238afSKris Buschelman if (size == 1) { 1128209238afSKris Buschelman ierr = MatSetType(A,MATSEQDENSE);CHKERRQ(ierr); 1129209238afSKris Buschelman } else { 1130209238afSKris Buschelman ierr = MatSetType(A,MATMPIDENSE);CHKERRQ(ierr); 1131209238afSKris Buschelman } 1132209238afSKris Buschelman PetscFunctionReturn(0); 1133209238afSKris Buschelman } 1134209238afSKris Buschelman EXTERN_C_END 1135209238afSKris Buschelman 11364a2ae208SSatish Balay #undef __FUNCT__ 11374a2ae208SSatish Balay #define __FUNCT__ "MatMPIDenseSetPreallocation" 1138273d9f13SBarry Smith /*@C 1139273d9f13SBarry Smith MatMPIDenseSetPreallocation - Sets the array used to store the matrix entries 1140273d9f13SBarry Smith 1141273d9f13SBarry Smith Not collective 1142273d9f13SBarry Smith 1143273d9f13SBarry Smith Input Parameters: 1144273d9f13SBarry Smith . A - the matrix 1145273d9f13SBarry Smith - data - optional location of matrix data. Set data=PETSC_NULL for PETSc 1146273d9f13SBarry Smith to control all matrix memory allocation. 1147273d9f13SBarry Smith 1148273d9f13SBarry Smith Notes: 1149273d9f13SBarry Smith The dense format is fully compatible with standard Fortran 77 1150273d9f13SBarry Smith storage by columns. 1151273d9f13SBarry Smith 1152273d9f13SBarry Smith The data input variable is intended primarily for Fortran programmers 1153273d9f13SBarry Smith who wish to allocate their own matrix memory space. Most users should 1154273d9f13SBarry Smith set data=PETSC_NULL. 1155273d9f13SBarry Smith 1156273d9f13SBarry Smith Level: intermediate 1157273d9f13SBarry Smith 1158273d9f13SBarry Smith .keywords: matrix,dense, parallel 1159273d9f13SBarry Smith 1160273d9f13SBarry Smith .seealso: MatCreate(), MatCreateSeqDense(), MatSetValues() 1161273d9f13SBarry Smith @*/ 1162be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatMPIDenseSetPreallocation(Mat mat,PetscScalar *data) 1163273d9f13SBarry Smith { 1164dfbe8321SBarry Smith PetscErrorCode ierr,(*f)(Mat,PetscScalar *); 1165273d9f13SBarry Smith 1166273d9f13SBarry Smith PetscFunctionBegin; 1167565567f0SKris Buschelman ierr = PetscObjectQueryFunction((PetscObject)mat,"MatMPIDenseSetPreallocation_C",(void (**)(void))&f);CHKERRQ(ierr); 1168a23d5eceSKris Buschelman if (f) { 1169a23d5eceSKris Buschelman ierr = (*f)(mat,data);CHKERRQ(ierr); 1170a23d5eceSKris Buschelman } 1171273d9f13SBarry Smith PetscFunctionReturn(0); 1172273d9f13SBarry Smith } 1173273d9f13SBarry Smith 11744a2ae208SSatish Balay #undef __FUNCT__ 11754a2ae208SSatish Balay #define __FUNCT__ "MatCreateMPIDense" 11768965ea79SLois Curfman McInnes /*@C 117739ddd567SLois Curfman McInnes MatCreateMPIDense - Creates a sparse parallel matrix in dense format. 11788965ea79SLois Curfman McInnes 1179db81eaa0SLois Curfman McInnes Collective on MPI_Comm 1180db81eaa0SLois Curfman McInnes 11818965ea79SLois Curfman McInnes Input Parameters: 1182db81eaa0SLois Curfman McInnes + comm - MPI communicator 11838965ea79SLois Curfman McInnes . m - number of local rows (or PETSC_DECIDE to have calculated if M is given) 1184db81eaa0SLois Curfman McInnes . n - number of local columns (or PETSC_DECIDE to have calculated if N is given) 11858965ea79SLois Curfman McInnes . M - number of global rows (or PETSC_DECIDE to have calculated if m is given) 1186db81eaa0SLois Curfman McInnes . N - number of global columns (or PETSC_DECIDE to have calculated if n is given) 11877f5ff6fdSBarry Smith - data - optional location of matrix data. Set data=PETSC_NULL (PETSC_NULL_SCALAR for Fortran users) for PETSc 1188dfc5480cSLois Curfman McInnes to control all matrix memory allocation. 11898965ea79SLois Curfman McInnes 11908965ea79SLois Curfman McInnes Output Parameter: 1191477f1c0bSLois Curfman McInnes . A - the matrix 11928965ea79SLois Curfman McInnes 1193b259b22eSLois Curfman McInnes Notes: 119439ddd567SLois Curfman McInnes The dense format is fully compatible with standard Fortran 77 119539ddd567SLois Curfman McInnes storage by columns. 11968965ea79SLois Curfman McInnes 119718f449edSLois Curfman McInnes The data input variable is intended primarily for Fortran programmers 119818f449edSLois Curfman McInnes who wish to allocate their own matrix memory space. Most users should 11997f5ff6fdSBarry Smith set data=PETSC_NULL (PETSC_NULL_SCALAR for Fortran users). 120018f449edSLois Curfman McInnes 12018965ea79SLois Curfman McInnes The user MUST specify either the local or global matrix dimensions 12028965ea79SLois Curfman McInnes (possibly both). 12038965ea79SLois Curfman McInnes 1204027ccd11SLois Curfman McInnes Level: intermediate 1205027ccd11SLois Curfman McInnes 120639ddd567SLois Curfman McInnes .keywords: matrix,dense, parallel 12078965ea79SLois Curfman McInnes 120839ddd567SLois Curfman McInnes .seealso: MatCreate(), MatCreateSeqDense(), MatSetValues() 12098965ea79SLois Curfman McInnes @*/ 1210be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIDense(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscScalar *data,Mat *A) 12118965ea79SLois Curfman McInnes { 12126849ba73SBarry Smith PetscErrorCode ierr; 121313f74950SBarry Smith PetscMPIInt size; 12148965ea79SLois Curfman McInnes 12153a40ed3dSBarry Smith PetscFunctionBegin; 1216f69a0ea3SMatthew Knepley ierr = MatCreate(comm,A);CHKERRQ(ierr); 1217f69a0ea3SMatthew Knepley ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr); 1218273d9f13SBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1219273d9f13SBarry Smith if (size > 1) { 1220273d9f13SBarry Smith ierr = MatSetType(*A,MATMPIDENSE);CHKERRQ(ierr); 1221273d9f13SBarry Smith ierr = MatMPIDenseSetPreallocation(*A,data);CHKERRQ(ierr); 1222273d9f13SBarry Smith } else { 1223273d9f13SBarry Smith ierr = MatSetType(*A,MATSEQDENSE);CHKERRQ(ierr); 1224273d9f13SBarry Smith ierr = MatSeqDenseSetPreallocation(*A,data);CHKERRQ(ierr); 12258c469469SLois Curfman McInnes } 12263a40ed3dSBarry Smith PetscFunctionReturn(0); 12278965ea79SLois Curfman McInnes } 12288965ea79SLois Curfman McInnes 12294a2ae208SSatish Balay #undef __FUNCT__ 12304a2ae208SSatish Balay #define __FUNCT__ "MatDuplicate_MPIDense" 12316849ba73SBarry Smith static PetscErrorCode MatDuplicate_MPIDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat) 12328965ea79SLois Curfman McInnes { 12338965ea79SLois Curfman McInnes Mat mat; 12343501a2bdSLois Curfman McInnes Mat_MPIDense *a,*oldmat = (Mat_MPIDense*)A->data; 1235dfbe8321SBarry Smith PetscErrorCode ierr; 12368965ea79SLois Curfman McInnes 12373a40ed3dSBarry Smith PetscFunctionBegin; 12388965ea79SLois Curfman McInnes *newmat = 0; 1239f69a0ea3SMatthew Knepley ierr = MatCreate(A->comm,&mat);CHKERRQ(ierr); 1240f69a0ea3SMatthew Knepley ierr = MatSetSizes(mat,A->m,A->n,A->M,A->N);CHKERRQ(ierr); 1241be5d1d56SKris Buschelman ierr = MatSetType(mat,A->type_name);CHKERRQ(ierr); 1242b0a32e0cSBarry Smith ierr = PetscNew(Mat_MPIDense,&a);CHKERRQ(ierr); 1243b0a32e0cSBarry Smith mat->data = (void*)a; 1244549d3d68SSatish Balay ierr = PetscMemcpy(mat->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 12453501a2bdSLois Curfman McInnes mat->factor = A->factor; 1246c456f294SBarry Smith mat->assembled = PETSC_TRUE; 1247273d9f13SBarry Smith mat->preallocated = PETSC_TRUE; 12488965ea79SLois Curfman McInnes 12498965ea79SLois Curfman McInnes a->rstart = oldmat->rstart; 12508965ea79SLois Curfman McInnes a->rend = oldmat->rend; 12518965ea79SLois Curfman McInnes a->size = oldmat->size; 12528965ea79SLois Curfman McInnes a->rank = oldmat->rank; 1253e0fa3b82SLois Curfman McInnes mat->insertmode = NOT_SET_VALUES; 1254b9b97703SBarry Smith a->nvec = oldmat->nvec; 12553782ba37SSatish Balay a->donotstash = oldmat->donotstash; 125613f74950SBarry Smith ierr = PetscMalloc((a->size+1)*sizeof(PetscInt),&a->rowners);CHKERRQ(ierr); 125752e6d16bSBarry Smith ierr = PetscLogObjectMemory(mat,(a->size+1)*sizeof(PetscInt)+sizeof(struct _p_Mat)+sizeof(Mat_MPIDense));CHKERRQ(ierr); 125813f74950SBarry Smith ierr = PetscMemcpy(a->rowners,oldmat->rowners,(a->size+1)*sizeof(PetscInt));CHKERRQ(ierr); 12598798bf22SSatish Balay ierr = MatStashCreate_Private(A->comm,1,&mat->stash);CHKERRQ(ierr); 12608965ea79SLois Curfman McInnes 1261329f5518SBarry Smith ierr = MatSetUpMultiply_MPIDense(mat);CHKERRQ(ierr); 12625609ef8eSBarry Smith ierr = MatDuplicate(oldmat->A,cpvalues,&a->A);CHKERRQ(ierr); 126352e6d16bSBarry Smith ierr = PetscLogObjectParent(mat,a->A);CHKERRQ(ierr); 12648965ea79SLois Curfman McInnes *newmat = mat; 12653a40ed3dSBarry Smith PetscFunctionReturn(0); 12668965ea79SLois Curfman McInnes } 12678965ea79SLois Curfman McInnes 1268e090d566SSatish Balay #include "petscsys.h" 12698965ea79SLois Curfman McInnes 12704a2ae208SSatish Balay #undef __FUNCT__ 12714a2ae208SSatish Balay #define __FUNCT__ "MatLoad_MPIDense_DenseInFile" 1272f69a0ea3SMatthew Knepley PetscErrorCode MatLoad_MPIDense_DenseInFile(MPI_Comm comm,PetscInt fd,PetscInt M,PetscInt N, MatType type,Mat *newmat) 127390ace30eSBarry Smith { 12746849ba73SBarry Smith PetscErrorCode ierr; 127532dcc486SBarry Smith PetscMPIInt rank,size; 127613f74950SBarry Smith PetscInt *rowners,i,m,nz,j; 127787828ca2SBarry Smith PetscScalar *array,*vals,*vals_ptr; 127890ace30eSBarry Smith MPI_Status status; 127990ace30eSBarry Smith 12803a40ed3dSBarry Smith PetscFunctionBegin; 1281d132466eSBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 1282d132466eSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 128390ace30eSBarry Smith 128490ace30eSBarry Smith /* determine ownership of all rows */ 128590ace30eSBarry Smith m = M/size + ((M % size) > rank); 128613f74950SBarry Smith ierr = PetscMalloc((size+2)*sizeof(PetscInt),&rowners);CHKERRQ(ierr); 128713f74950SBarry Smith ierr = MPI_Allgather(&m,1,MPIU_INT,rowners+1,1,MPIU_INT,comm);CHKERRQ(ierr); 128890ace30eSBarry Smith rowners[0] = 0; 128990ace30eSBarry Smith for (i=2; i<=size; i++) { 129090ace30eSBarry Smith rowners[i] += rowners[i-1]; 129190ace30eSBarry Smith } 129290ace30eSBarry Smith 1293f69a0ea3SMatthew Knepley ierr = MatCreate(comm,newmat);CHKERRQ(ierr); 1294f69a0ea3SMatthew Knepley ierr = MatSetSizes(*newmat,m,PETSC_DECIDE,M,N);CHKERRQ(ierr); 1295be5d1d56SKris Buschelman ierr = MatSetType(*newmat,type);CHKERRQ(ierr); 1296878740d9SKris Buschelman ierr = MatMPIDenseSetPreallocation(*newmat,PETSC_NULL);CHKERRQ(ierr); 129790ace30eSBarry Smith ierr = MatGetArray(*newmat,&array);CHKERRQ(ierr); 129890ace30eSBarry Smith 129990ace30eSBarry Smith if (!rank) { 130087828ca2SBarry Smith ierr = PetscMalloc(m*N*sizeof(PetscScalar),&vals);CHKERRQ(ierr); 130190ace30eSBarry Smith 130290ace30eSBarry Smith /* read in my part of the matrix numerical values */ 13030752156aSBarry Smith ierr = PetscBinaryRead(fd,vals,m*N,PETSC_SCALAR);CHKERRQ(ierr); 130490ace30eSBarry Smith 130590ace30eSBarry Smith /* insert into matrix-by row (this is why cannot directly read into array */ 130690ace30eSBarry Smith vals_ptr = vals; 130790ace30eSBarry Smith for (i=0; i<m; i++) { 130890ace30eSBarry Smith for (j=0; j<N; j++) { 130990ace30eSBarry Smith array[i + j*m] = *vals_ptr++; 131090ace30eSBarry Smith } 131190ace30eSBarry Smith } 131290ace30eSBarry Smith 131390ace30eSBarry Smith /* read in other processors and ship out */ 131490ace30eSBarry Smith for (i=1; i<size; i++) { 131590ace30eSBarry Smith nz = (rowners[i+1] - rowners[i])*N; 13160752156aSBarry Smith ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr); 1317ca161407SBarry Smith ierr = MPI_Send(vals,nz,MPIU_SCALAR,i,(*newmat)->tag,comm);CHKERRQ(ierr); 131890ace30eSBarry Smith } 13193a40ed3dSBarry Smith } else { 132090ace30eSBarry Smith /* receive numeric values */ 132187828ca2SBarry Smith ierr = PetscMalloc(m*N*sizeof(PetscScalar),&vals);CHKERRQ(ierr); 132290ace30eSBarry Smith 132390ace30eSBarry Smith /* receive message of values*/ 1324ca161407SBarry Smith ierr = MPI_Recv(vals,m*N,MPIU_SCALAR,0,(*newmat)->tag,comm,&status);CHKERRQ(ierr); 132590ace30eSBarry Smith 132690ace30eSBarry Smith /* insert into matrix-by row (this is why cannot directly read into array */ 132790ace30eSBarry Smith vals_ptr = vals; 132890ace30eSBarry Smith for (i=0; i<m; i++) { 132990ace30eSBarry Smith for (j=0; j<N; j++) { 133090ace30eSBarry Smith array[i + j*m] = *vals_ptr++; 133190ace30eSBarry Smith } 133290ace30eSBarry Smith } 133390ace30eSBarry Smith } 1334606d414cSSatish Balay ierr = PetscFree(rowners);CHKERRQ(ierr); 1335606d414cSSatish Balay ierr = PetscFree(vals);CHKERRQ(ierr); 13366d4a8577SBarry Smith ierr = MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 13376d4a8577SBarry Smith ierr = MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 13383a40ed3dSBarry Smith PetscFunctionReturn(0); 133990ace30eSBarry Smith } 134090ace30eSBarry Smith 13414a2ae208SSatish Balay #undef __FUNCT__ 13424a2ae208SSatish Balay #define __FUNCT__ "MatLoad_MPIDense" 1343f69a0ea3SMatthew Knepley PetscErrorCode MatLoad_MPIDense(PetscViewer viewer, MatType type,Mat *newmat) 13448965ea79SLois Curfman McInnes { 13458965ea79SLois Curfman McInnes Mat A; 134687828ca2SBarry Smith PetscScalar *vals,*svals; 134719bcc07fSBarry Smith MPI_Comm comm = ((PetscObject)viewer)->comm; 13488965ea79SLois Curfman McInnes MPI_Status status; 134913f74950SBarry Smith PetscMPIInt rank,size,tag = ((PetscObject)viewer)->tag,*rowners,*sndcounts,m,maxnz; 135013f74950SBarry Smith PetscInt header[4],*rowlengths = 0,M,N,*cols; 135113f74950SBarry Smith PetscInt *ourlens,*procsnz = 0,*offlens,jj,*mycols,*smycols; 135213f74950SBarry Smith PetscInt i,nz,j,rstart,rend; 135313f74950SBarry Smith int fd; 13546849ba73SBarry Smith PetscErrorCode ierr; 13558965ea79SLois Curfman McInnes 13563a40ed3dSBarry Smith PetscFunctionBegin; 1357d132466eSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1358d132466eSBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 13598965ea79SLois Curfman McInnes if (!rank) { 1360b0a32e0cSBarry Smith ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 13610752156aSBarry Smith ierr = PetscBinaryRead(fd,(char *)header,4,PETSC_INT);CHKERRQ(ierr); 1362552e946dSBarry Smith if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"not matrix object"); 13638965ea79SLois Curfman McInnes } 13648965ea79SLois Curfman McInnes 136513f74950SBarry Smith ierr = MPI_Bcast(header+1,3,MPIU_INT,0,comm);CHKERRQ(ierr); 136690ace30eSBarry Smith M = header[1]; N = header[2]; nz = header[3]; 136790ace30eSBarry Smith 136890ace30eSBarry Smith /* 136990ace30eSBarry Smith Handle case where matrix is stored on disk as a dense matrix 137090ace30eSBarry Smith */ 137190ace30eSBarry Smith if (nz == MATRIX_BINARY_FORMAT_DENSE) { 1372be5d1d56SKris Buschelman ierr = MatLoad_MPIDense_DenseInFile(comm,fd,M,N,type,newmat);CHKERRQ(ierr); 13733a40ed3dSBarry Smith PetscFunctionReturn(0); 137490ace30eSBarry Smith } 137590ace30eSBarry Smith 13768965ea79SLois Curfman McInnes /* determine ownership of all rows */ 13778965ea79SLois Curfman McInnes m = M/size + ((M % size) > rank); 137813f74950SBarry Smith ierr = PetscMalloc((size+2)*sizeof(PetscInt),&rowners);CHKERRQ(ierr); 1379ca161407SBarry Smith ierr = MPI_Allgather(&m,1,MPI_INT,rowners+1,1,MPI_INT,comm);CHKERRQ(ierr); 13808965ea79SLois Curfman McInnes rowners[0] = 0; 13818965ea79SLois Curfman McInnes for (i=2; i<=size; i++) { 13828965ea79SLois Curfman McInnes rowners[i] += rowners[i-1]; 13838965ea79SLois Curfman McInnes } 13848965ea79SLois Curfman McInnes rstart = rowners[rank]; 13858965ea79SLois Curfman McInnes rend = rowners[rank+1]; 13868965ea79SLois Curfman McInnes 13878965ea79SLois Curfman McInnes /* distribute row lengths to all processors */ 138813f74950SBarry Smith ierr = PetscMalloc(2*(rend-rstart+1)*sizeof(PetscInt),&ourlens);CHKERRQ(ierr); 13898965ea79SLois Curfman McInnes offlens = ourlens + (rend-rstart); 13908965ea79SLois Curfman McInnes if (!rank) { 139113f74950SBarry Smith ierr = PetscMalloc(M*sizeof(PetscInt),&rowlengths);CHKERRQ(ierr); 13920752156aSBarry Smith ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr); 139313f74950SBarry Smith ierr = PetscMalloc(size*sizeof(PetscMPIInt),&sndcounts);CHKERRQ(ierr); 13948965ea79SLois Curfman McInnes for (i=0; i<size; i++) sndcounts[i] = rowners[i+1] - rowners[i]; 13951466f1e1SBarry Smith ierr = MPI_Scatterv(rowlengths,sndcounts,rowners,MPIU_INT,ourlens,rend-rstart,MPIU_INT,0,comm);CHKERRQ(ierr); 1396606d414cSSatish Balay ierr = PetscFree(sndcounts);CHKERRQ(ierr); 1397ca161407SBarry Smith } else { 13981466f1e1SBarry Smith ierr = MPI_Scatterv(0,0,0,MPIU_INT,ourlens,rend-rstart,MPIU_INT,0,comm);CHKERRQ(ierr); 13998965ea79SLois Curfman McInnes } 14008965ea79SLois Curfman McInnes 14018965ea79SLois Curfman McInnes if (!rank) { 14028965ea79SLois Curfman McInnes /* calculate the number of nonzeros on each processor */ 140313f74950SBarry Smith ierr = PetscMalloc(size*sizeof(PetscInt),&procsnz);CHKERRQ(ierr); 140413f74950SBarry Smith ierr = PetscMemzero(procsnz,size*sizeof(PetscInt));CHKERRQ(ierr); 14058965ea79SLois Curfman McInnes for (i=0; i<size; i++) { 14068965ea79SLois Curfman McInnes for (j=rowners[i]; j< rowners[i+1]; j++) { 14078965ea79SLois Curfman McInnes procsnz[i] += rowlengths[j]; 14088965ea79SLois Curfman McInnes } 14098965ea79SLois Curfman McInnes } 1410606d414cSSatish Balay ierr = PetscFree(rowlengths);CHKERRQ(ierr); 14118965ea79SLois Curfman McInnes 14128965ea79SLois Curfman McInnes /* determine max buffer needed and allocate it */ 14138965ea79SLois Curfman McInnes maxnz = 0; 14148965ea79SLois Curfman McInnes for (i=0; i<size; i++) { 14150452661fSBarry Smith maxnz = PetscMax(maxnz,procsnz[i]); 14168965ea79SLois Curfman McInnes } 141713f74950SBarry Smith ierr = PetscMalloc(maxnz*sizeof(PetscInt),&cols);CHKERRQ(ierr); 14188965ea79SLois Curfman McInnes 14198965ea79SLois Curfman McInnes /* read in my part of the matrix column indices */ 14208965ea79SLois Curfman McInnes nz = procsnz[0]; 142113f74950SBarry Smith ierr = PetscMalloc(nz*sizeof(PetscInt),&mycols);CHKERRQ(ierr); 14220752156aSBarry Smith ierr = PetscBinaryRead(fd,mycols,nz,PETSC_INT);CHKERRQ(ierr); 14238965ea79SLois Curfman McInnes 14248965ea79SLois Curfman McInnes /* read in every one elses and ship off */ 14258965ea79SLois Curfman McInnes for (i=1; i<size; i++) { 14268965ea79SLois Curfman McInnes nz = procsnz[i]; 14270752156aSBarry Smith ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT);CHKERRQ(ierr); 142813f74950SBarry Smith ierr = MPI_Send(cols,nz,MPIU_INT,i,tag,comm);CHKERRQ(ierr); 14298965ea79SLois Curfman McInnes } 1430606d414cSSatish Balay ierr = PetscFree(cols);CHKERRQ(ierr); 14313a40ed3dSBarry Smith } else { 14328965ea79SLois Curfman McInnes /* determine buffer space needed for message */ 14338965ea79SLois Curfman McInnes nz = 0; 14348965ea79SLois Curfman McInnes for (i=0; i<m; i++) { 14358965ea79SLois Curfman McInnes nz += ourlens[i]; 14368965ea79SLois Curfman McInnes } 143713f74950SBarry Smith ierr = PetscMalloc((nz+1)*sizeof(PetscInt),&mycols);CHKERRQ(ierr); 14388965ea79SLois Curfman McInnes 14398965ea79SLois Curfman McInnes /* receive message of column indices*/ 144013f74950SBarry Smith ierr = MPI_Recv(mycols,nz,MPIU_INT,0,tag,comm,&status);CHKERRQ(ierr); 144113f74950SBarry Smith ierr = MPI_Get_count(&status,MPIU_INT,&maxnz);CHKERRQ(ierr); 144229bbc08cSBarry Smith if (maxnz != nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"something is wrong with file"); 14438965ea79SLois Curfman McInnes } 14448965ea79SLois Curfman McInnes 14458965ea79SLois Curfman McInnes /* loop over local rows, determining number of off diagonal entries */ 144613f74950SBarry Smith ierr = PetscMemzero(offlens,m*sizeof(PetscInt));CHKERRQ(ierr); 14478965ea79SLois Curfman McInnes jj = 0; 14488965ea79SLois Curfman McInnes for (i=0; i<m; i++) { 14498965ea79SLois Curfman McInnes for (j=0; j<ourlens[i]; j++) { 14508965ea79SLois Curfman McInnes if (mycols[jj] < rstart || mycols[jj] >= rend) offlens[i]++; 14518965ea79SLois Curfman McInnes jj++; 14528965ea79SLois Curfman McInnes } 14538965ea79SLois Curfman McInnes } 14548965ea79SLois Curfman McInnes 14558965ea79SLois Curfman McInnes /* create our matrix */ 14568965ea79SLois Curfman McInnes for (i=0; i<m; i++) { 14578965ea79SLois Curfman McInnes ourlens[i] -= offlens[i]; 14588965ea79SLois Curfman McInnes } 1459f69a0ea3SMatthew Knepley ierr = MatCreate(comm,newmat);CHKERRQ(ierr); 1460f69a0ea3SMatthew Knepley ierr = MatSetSizes(*newmat,m,PETSC_DECIDE,M,N);CHKERRQ(ierr); 1461be5d1d56SKris Buschelman ierr = MatSetType(*newmat,type);CHKERRQ(ierr); 1462878740d9SKris Buschelman ierr = MatMPIDenseSetPreallocation(*newmat,PETSC_NULL);CHKERRQ(ierr); 14638965ea79SLois Curfman McInnes A = *newmat; 14648965ea79SLois Curfman McInnes for (i=0; i<m; i++) { 14658965ea79SLois Curfman McInnes ourlens[i] += offlens[i]; 14668965ea79SLois Curfman McInnes } 14678965ea79SLois Curfman McInnes 14688965ea79SLois Curfman McInnes if (!rank) { 146987828ca2SBarry Smith ierr = PetscMalloc(maxnz*sizeof(PetscScalar),&vals);CHKERRQ(ierr); 14708965ea79SLois Curfman McInnes 14718965ea79SLois Curfman McInnes /* read in my part of the matrix numerical values */ 14728965ea79SLois Curfman McInnes nz = procsnz[0]; 14730752156aSBarry Smith ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr); 14748965ea79SLois Curfman McInnes 14758965ea79SLois Curfman McInnes /* insert into matrix */ 14768965ea79SLois Curfman McInnes jj = rstart; 14778965ea79SLois Curfman McInnes smycols = mycols; 14788965ea79SLois Curfman McInnes svals = vals; 14798965ea79SLois Curfman McInnes for (i=0; i<m; i++) { 14808965ea79SLois Curfman McInnes ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr); 14818965ea79SLois Curfman McInnes smycols += ourlens[i]; 14828965ea79SLois Curfman McInnes svals += ourlens[i]; 14838965ea79SLois Curfman McInnes jj++; 14848965ea79SLois Curfman McInnes } 14858965ea79SLois Curfman McInnes 14868965ea79SLois Curfman McInnes /* read in other processors and ship out */ 14878965ea79SLois Curfman McInnes for (i=1; i<size; i++) { 14888965ea79SLois Curfman McInnes nz = procsnz[i]; 14890752156aSBarry Smith ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr); 1490ca161407SBarry Smith ierr = MPI_Send(vals,nz,MPIU_SCALAR,i,A->tag,comm);CHKERRQ(ierr); 14918965ea79SLois Curfman McInnes } 1492606d414cSSatish Balay ierr = PetscFree(procsnz);CHKERRQ(ierr); 14933a40ed3dSBarry Smith } else { 14948965ea79SLois Curfman McInnes /* receive numeric values */ 149584d528b1SBarry Smith ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&vals);CHKERRQ(ierr); 14968965ea79SLois Curfman McInnes 14978965ea79SLois Curfman McInnes /* receive message of values*/ 1498ca161407SBarry Smith ierr = MPI_Recv(vals,nz,MPIU_SCALAR,0,A->tag,comm,&status);CHKERRQ(ierr); 1499ca161407SBarry Smith ierr = MPI_Get_count(&status,MPIU_SCALAR,&maxnz);CHKERRQ(ierr); 150029bbc08cSBarry Smith if (maxnz != nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"something is wrong with file"); 15018965ea79SLois Curfman McInnes 15028965ea79SLois Curfman McInnes /* insert into matrix */ 15038965ea79SLois Curfman McInnes jj = rstart; 15048965ea79SLois Curfman McInnes smycols = mycols; 15058965ea79SLois Curfman McInnes svals = vals; 15068965ea79SLois Curfman McInnes for (i=0; i<m; i++) { 15078965ea79SLois Curfman McInnes ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr); 15088965ea79SLois Curfman McInnes smycols += ourlens[i]; 15098965ea79SLois Curfman McInnes svals += ourlens[i]; 15108965ea79SLois Curfman McInnes jj++; 15118965ea79SLois Curfman McInnes } 15128965ea79SLois Curfman McInnes } 1513606d414cSSatish Balay ierr = PetscFree(ourlens);CHKERRQ(ierr); 1514606d414cSSatish Balay ierr = PetscFree(vals);CHKERRQ(ierr); 1515606d414cSSatish Balay ierr = PetscFree(mycols);CHKERRQ(ierr); 1516606d414cSSatish Balay ierr = PetscFree(rowners);CHKERRQ(ierr); 15178965ea79SLois Curfman McInnes 15186d4a8577SBarry Smith ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 15196d4a8577SBarry Smith ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 15203a40ed3dSBarry Smith PetscFunctionReturn(0); 15218965ea79SLois Curfman McInnes } 152290ace30eSBarry Smith 152390ace30eSBarry Smith 152490ace30eSBarry Smith 152590ace30eSBarry Smith 152690ace30eSBarry Smith 1527