1be1d678aSKris Buschelman 2ed3cc1f0SBarry Smith /* 3ed3cc1f0SBarry Smith Basic functions for basic parallel dense matrices. 4ed3cc1f0SBarry Smith */ 5ed3cc1f0SBarry Smith 604fea9ffSBarry Smith 7c6db04a5SJed Brown #include <../src/mat/impls/dense/mpi/mpidense.h> /*I "petscmat.h" I*/ 88949adfdSHong Zhang #include <../src/mat/impls/aij/mpi/mpiaij.h> 9baa3c1c6SHong Zhang #include <petscblaslapack.h> 108965ea79SLois Curfman McInnes 11ab92ecdeSBarry Smith /*@ 12ab92ecdeSBarry Smith 13ab92ecdeSBarry Smith MatDenseGetLocalMatrix - For a MATMPIDENSE or MATSEQDENSE matrix returns the sequential 14ab92ecdeSBarry Smith matrix that represents the operator. For sequential matrices it returns itself. 15ab92ecdeSBarry Smith 16ab92ecdeSBarry Smith Input Parameter: 17ab92ecdeSBarry Smith . A - the Seq or MPI dense matrix 18ab92ecdeSBarry Smith 19ab92ecdeSBarry Smith Output Parameter: 20ab92ecdeSBarry Smith . B - the inner matrix 21ab92ecdeSBarry Smith 228e6c10adSSatish Balay Level: intermediate 238e6c10adSSatish Balay 24ab92ecdeSBarry Smith @*/ 25ab92ecdeSBarry Smith PetscErrorCode MatDenseGetLocalMatrix(Mat A,Mat *B) 26ab92ecdeSBarry Smith { 27ab92ecdeSBarry Smith Mat_MPIDense *mat = (Mat_MPIDense*)A->data; 28ab92ecdeSBarry Smith PetscErrorCode ierr; 29ace3abfcSBarry Smith PetscBool flg; 30ab92ecdeSBarry Smith 31ab92ecdeSBarry Smith PetscFunctionBegin; 32251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIDENSE,&flg);CHKERRQ(ierr); 332205254eSKarl Rupp if (flg) *B = mat->A; 342205254eSKarl Rupp else *B = A; 35ab92ecdeSBarry Smith PetscFunctionReturn(0); 36ab92ecdeSBarry Smith } 37ab92ecdeSBarry Smith 38ba8c8a56SBarry Smith PetscErrorCode MatGetRow_MPIDense(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v) 39ba8c8a56SBarry Smith { 40ba8c8a56SBarry Smith Mat_MPIDense *mat = (Mat_MPIDense*)A->data; 41ba8c8a56SBarry Smith PetscErrorCode ierr; 42d0f46423SBarry Smith PetscInt lrow,rstart = A->rmap->rstart,rend = A->rmap->rend; 43ba8c8a56SBarry Smith 44ba8c8a56SBarry Smith PetscFunctionBegin; 45e7e72b3dSBarry Smith if (row < rstart || row >= rend) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"only local rows"); 46ba8c8a56SBarry Smith lrow = row - rstart; 47ba8c8a56SBarry Smith ierr = MatGetRow(mat->A,lrow,nz,(const PetscInt**)idx,(const PetscScalar**)v);CHKERRQ(ierr); 48ba8c8a56SBarry Smith PetscFunctionReturn(0); 49ba8c8a56SBarry Smith } 50ba8c8a56SBarry Smith 51ba8c8a56SBarry Smith PetscErrorCode MatRestoreRow_MPIDense(Mat mat,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v) 52ba8c8a56SBarry Smith { 53ba8c8a56SBarry Smith PetscErrorCode ierr; 54ba8c8a56SBarry Smith 55ba8c8a56SBarry Smith PetscFunctionBegin; 56ba8c8a56SBarry Smith if (idx) {ierr = PetscFree(*idx);CHKERRQ(ierr);} 57ba8c8a56SBarry Smith if (v) {ierr = PetscFree(*v);CHKERRQ(ierr);} 58ba8c8a56SBarry Smith PetscFunctionReturn(0); 59ba8c8a56SBarry Smith } 60ba8c8a56SBarry Smith 6111bd1e4dSLisandro Dalcin PetscErrorCode MatGetDiagonalBlock_MPIDense(Mat A,Mat *a) 620de54da6SSatish Balay { 630de54da6SSatish Balay Mat_MPIDense *mdn = (Mat_MPIDense*)A->data; 646849ba73SBarry Smith PetscErrorCode ierr; 65d0f46423SBarry Smith PetscInt m = A->rmap->n,rstart = A->rmap->rstart; 6687828ca2SBarry Smith PetscScalar *array; 670de54da6SSatish Balay MPI_Comm comm; 6811bd1e4dSLisandro Dalcin Mat B; 690de54da6SSatish Balay 700de54da6SSatish Balay PetscFunctionBegin; 71e32f2f54SBarry Smith if (A->rmap->N != A->cmap->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only square matrices supported."); 720de54da6SSatish Balay 7311bd1e4dSLisandro Dalcin ierr = PetscObjectQuery((PetscObject)A,"DiagonalBlock",(PetscObject*)&B);CHKERRQ(ierr); 7411bd1e4dSLisandro Dalcin if (!B) { 750de54da6SSatish Balay ierr = PetscObjectGetComm((PetscObject)(mdn->A),&comm);CHKERRQ(ierr); 7611bd1e4dSLisandro Dalcin ierr = MatCreate(comm,&B);CHKERRQ(ierr); 7711bd1e4dSLisandro Dalcin ierr = MatSetSizes(B,m,m,m,m);CHKERRQ(ierr); 7811bd1e4dSLisandro Dalcin ierr = MatSetType(B,((PetscObject)mdn->A)->type_name);CHKERRQ(ierr); 798c778c55SBarry Smith ierr = MatDenseGetArray(mdn->A,&array);CHKERRQ(ierr); 8011bd1e4dSLisandro Dalcin ierr = MatSeqDenseSetPreallocation(B,array+m*rstart);CHKERRQ(ierr); 818c778c55SBarry Smith ierr = MatDenseRestoreArray(mdn->A,&array);CHKERRQ(ierr); 8211bd1e4dSLisandro Dalcin ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 8311bd1e4dSLisandro Dalcin ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 8411bd1e4dSLisandro Dalcin ierr = PetscObjectCompose((PetscObject)A,"DiagonalBlock",(PetscObject)B);CHKERRQ(ierr); 8511bd1e4dSLisandro Dalcin *a = B; 8611bd1e4dSLisandro Dalcin ierr = MatDestroy(&B);CHKERRQ(ierr); 872205254eSKarl Rupp } else *a = B; 880de54da6SSatish Balay PetscFunctionReturn(0); 890de54da6SSatish Balay } 900de54da6SSatish Balay 9113f74950SBarry Smith PetscErrorCode MatSetValues_MPIDense(Mat mat,PetscInt m,const PetscInt idxm[],PetscInt n,const PetscInt idxn[],const PetscScalar v[],InsertMode addv) 928965ea79SLois Curfman McInnes { 9339b7565bSBarry Smith Mat_MPIDense *A = (Mat_MPIDense*)mat->data; 94dfbe8321SBarry Smith PetscErrorCode ierr; 95d0f46423SBarry Smith PetscInt i,j,rstart = mat->rmap->rstart,rend = mat->rmap->rend,row; 96ace3abfcSBarry Smith PetscBool roworiented = A->roworiented; 978965ea79SLois Curfman McInnes 983a40ed3dSBarry Smith PetscFunctionBegin; 998965ea79SLois Curfman McInnes for (i=0; i<m; i++) { 1005ef9f2a5SBarry Smith if (idxm[i] < 0) continue; 101e32f2f54SBarry Smith if (idxm[i] >= mat->rmap->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row too large"); 1028965ea79SLois Curfman McInnes if (idxm[i] >= rstart && idxm[i] < rend) { 1038965ea79SLois Curfman McInnes row = idxm[i] - rstart; 10439b7565bSBarry Smith if (roworiented) { 10539b7565bSBarry Smith ierr = MatSetValues(A->A,1,&row,n,idxn,v+i*n,addv);CHKERRQ(ierr); 1063a40ed3dSBarry Smith } else { 1078965ea79SLois Curfman McInnes for (j=0; j<n; j++) { 1085ef9f2a5SBarry Smith if (idxn[j] < 0) continue; 109e32f2f54SBarry Smith if (idxn[j] >= mat->cmap->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column too large"); 11039b7565bSBarry Smith ierr = MatSetValues(A->A,1,&row,1,&idxn[j],v+i+j*m,addv);CHKERRQ(ierr); 11139b7565bSBarry Smith } 1128965ea79SLois Curfman McInnes } 1132205254eSKarl Rupp } else if (!A->donotstash) { 1145080c13bSMatthew G Knepley mat->assembled = PETSC_FALSE; 11539b7565bSBarry Smith if (roworiented) { 116b400d20cSBarry Smith ierr = MatStashValuesRow_Private(&mat->stash,idxm[i],n,idxn,v+i*n,PETSC_FALSE);CHKERRQ(ierr); 117d36fbae8SSatish Balay } else { 118b400d20cSBarry Smith ierr = MatStashValuesCol_Private(&mat->stash,idxm[i],n,idxn,v+i,m,PETSC_FALSE);CHKERRQ(ierr); 11939b7565bSBarry Smith } 120b49de8d1SLois Curfman McInnes } 121b49de8d1SLois Curfman McInnes } 1223a40ed3dSBarry Smith PetscFunctionReturn(0); 123b49de8d1SLois Curfman McInnes } 124b49de8d1SLois Curfman McInnes 12513f74950SBarry Smith PetscErrorCode MatGetValues_MPIDense(Mat mat,PetscInt m,const PetscInt idxm[],PetscInt n,const PetscInt idxn[],PetscScalar v[]) 126b49de8d1SLois Curfman McInnes { 127b49de8d1SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 128dfbe8321SBarry Smith PetscErrorCode ierr; 129d0f46423SBarry Smith PetscInt i,j,rstart = mat->rmap->rstart,rend = mat->rmap->rend,row; 130b49de8d1SLois Curfman McInnes 1313a40ed3dSBarry Smith PetscFunctionBegin; 132b49de8d1SLois Curfman McInnes for (i=0; i<m; i++) { 133e32f2f54SBarry Smith if (idxm[i] < 0) continue; /* SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative row"); */ 134e32f2f54SBarry Smith if (idxm[i] >= mat->rmap->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row too large"); 135b49de8d1SLois Curfman McInnes if (idxm[i] >= rstart && idxm[i] < rend) { 136b49de8d1SLois Curfman McInnes row = idxm[i] - rstart; 137b49de8d1SLois Curfman McInnes for (j=0; j<n; j++) { 138e32f2f54SBarry Smith if (idxn[j] < 0) continue; /* SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative column"); */ 139e7e72b3dSBarry Smith if (idxn[j] >= mat->cmap->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column too large"); 140b49de8d1SLois Curfman McInnes ierr = MatGetValues(mdn->A,1,&row,1,&idxn[j],v+i*n+j);CHKERRQ(ierr); 141b49de8d1SLois Curfman McInnes } 142e7e72b3dSBarry Smith } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only local values currently supported"); 1438965ea79SLois Curfman McInnes } 1443a40ed3dSBarry Smith PetscFunctionReturn(0); 1458965ea79SLois Curfman McInnes } 1468965ea79SLois Curfman McInnes 147d3042a70SBarry Smith static PetscErrorCode MatDenseGetArray_MPIDense(Mat A,PetscScalar *array[]) 148ff14e315SSatish Balay { 149ff14e315SSatish Balay Mat_MPIDense *a = (Mat_MPIDense*)A->data; 150dfbe8321SBarry Smith PetscErrorCode ierr; 151ff14e315SSatish Balay 1523a40ed3dSBarry Smith PetscFunctionBegin; 1538c778c55SBarry Smith ierr = MatDenseGetArray(a->A,array);CHKERRQ(ierr); 1543a40ed3dSBarry Smith PetscFunctionReturn(0); 155ff14e315SSatish Balay } 156ff14e315SSatish Balay 157d3042a70SBarry Smith static PetscErrorCode MatDensePlaceArray_MPIDense(Mat A,const PetscScalar array[]) 158d3042a70SBarry Smith { 159d3042a70SBarry Smith Mat_MPIDense *a = (Mat_MPIDense*)A->data; 160d3042a70SBarry Smith PetscErrorCode ierr; 161d3042a70SBarry Smith 162d3042a70SBarry Smith PetscFunctionBegin; 163d3042a70SBarry Smith ierr = MatDensePlaceArray(a->A,array);CHKERRQ(ierr); 164d3042a70SBarry Smith PetscFunctionReturn(0); 165d3042a70SBarry Smith } 166d3042a70SBarry Smith 167d3042a70SBarry Smith static PetscErrorCode MatDenseResetArray_MPIDense(Mat A) 168d3042a70SBarry Smith { 169d3042a70SBarry Smith Mat_MPIDense *a = (Mat_MPIDense*)A->data; 170d3042a70SBarry Smith PetscErrorCode ierr; 171d3042a70SBarry Smith 172d3042a70SBarry Smith PetscFunctionBegin; 173d3042a70SBarry Smith ierr = MatDenseResetArray(a->A);CHKERRQ(ierr); 174d3042a70SBarry Smith PetscFunctionReturn(0); 175d3042a70SBarry Smith } 176d3042a70SBarry Smith 1777dae84e0SHong Zhang static PetscErrorCode MatCreateSubMatrix_MPIDense(Mat A,IS isrow,IS iscol,MatReuse scall,Mat *B) 178ca3fa75bSLois Curfman McInnes { 179ca3fa75bSLois Curfman McInnes Mat_MPIDense *mat = (Mat_MPIDense*)A->data,*newmatd; 180ca3fa75bSLois Curfman McInnes Mat_SeqDense *lmat = (Mat_SeqDense*)mat->A->data; 1816849ba73SBarry Smith PetscErrorCode ierr; 1824aa3045dSJed Brown PetscInt i,j,rstart,rend,nrows,ncols,Ncols,nlrows,nlcols; 1835d0c19d7SBarry Smith const PetscInt *irow,*icol; 18487828ca2SBarry Smith PetscScalar *av,*bv,*v = lmat->v; 185ca3fa75bSLois Curfman McInnes Mat newmat; 1864aa3045dSJed Brown IS iscol_local; 187ca3fa75bSLois Curfman McInnes 188ca3fa75bSLois Curfman McInnes PetscFunctionBegin; 1894aa3045dSJed Brown ierr = ISAllGather(iscol,&iscol_local);CHKERRQ(ierr); 190ca3fa75bSLois Curfman McInnes ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr); 1914aa3045dSJed Brown ierr = ISGetIndices(iscol_local,&icol);CHKERRQ(ierr); 192b9b97703SBarry Smith ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr); 193b9b97703SBarry Smith ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr); 1944aa3045dSJed Brown ierr = ISGetSize(iscol,&Ncols);CHKERRQ(ierr); /* global number of columns, size of iscol_local */ 195ca3fa75bSLois Curfman McInnes 196ca3fa75bSLois Curfman McInnes /* No parallel redistribution currently supported! Should really check each index set 1977eba5e9cSLois Curfman McInnes to comfirm that it is OK. ... Currently supports only submatrix same partitioning as 1987eba5e9cSLois Curfman McInnes original matrix! */ 199ca3fa75bSLois Curfman McInnes 200ca3fa75bSLois Curfman McInnes ierr = MatGetLocalSize(A,&nlrows,&nlcols);CHKERRQ(ierr); 2017eba5e9cSLois Curfman McInnes ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr); 202ca3fa75bSLois Curfman McInnes 203ca3fa75bSLois Curfman McInnes /* Check submatrix call */ 204ca3fa75bSLois Curfman McInnes if (scall == MAT_REUSE_MATRIX) { 205e32f2f54SBarry Smith /* SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Reused submatrix wrong size"); */ 2067eba5e9cSLois Curfman McInnes /* Really need to test rows and column sizes! */ 207ca3fa75bSLois Curfman McInnes newmat = *B; 208ca3fa75bSLois Curfman McInnes } else { 209ca3fa75bSLois Curfman McInnes /* Create and fill new matrix */ 210ce94432eSBarry Smith ierr = MatCreate(PetscObjectComm((PetscObject)A),&newmat);CHKERRQ(ierr); 2114aa3045dSJed Brown ierr = MatSetSizes(newmat,nrows,ncols,PETSC_DECIDE,Ncols);CHKERRQ(ierr); 2127adad957SLisandro Dalcin ierr = MatSetType(newmat,((PetscObject)A)->type_name);CHKERRQ(ierr); 2130298fd71SBarry Smith ierr = MatMPIDenseSetPreallocation(newmat,NULL);CHKERRQ(ierr); 214ca3fa75bSLois Curfman McInnes } 215ca3fa75bSLois Curfman McInnes 216ca3fa75bSLois Curfman McInnes /* Now extract the data pointers and do the copy, column at a time */ 217ca3fa75bSLois Curfman McInnes newmatd = (Mat_MPIDense*)newmat->data; 218ca3fa75bSLois Curfman McInnes bv = ((Mat_SeqDense*)newmatd->A->data)->v; 219ca3fa75bSLois Curfman McInnes 2204aa3045dSJed Brown for (i=0; i<Ncols; i++) { 22125a33276SHong Zhang av = v + ((Mat_SeqDense*)mat->A->data)->lda*icol[i]; 222ca3fa75bSLois Curfman McInnes for (j=0; j<nrows; j++) { 2237eba5e9cSLois Curfman McInnes *bv++ = av[irow[j] - rstart]; 224ca3fa75bSLois Curfman McInnes } 225ca3fa75bSLois Curfman McInnes } 226ca3fa75bSLois Curfman McInnes 227ca3fa75bSLois Curfman McInnes /* Assemble the matrices so that the correct flags are set */ 228ca3fa75bSLois Curfman McInnes ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 229ca3fa75bSLois Curfman McInnes ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 230ca3fa75bSLois Curfman McInnes 231ca3fa75bSLois Curfman McInnes /* Free work space */ 232ca3fa75bSLois Curfman McInnes ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr); 2335bdf786aSShri Abhyankar ierr = ISRestoreIndices(iscol_local,&icol);CHKERRQ(ierr); 23432bb1f2dSLisandro Dalcin ierr = ISDestroy(&iscol_local);CHKERRQ(ierr); 235ca3fa75bSLois Curfman McInnes *B = newmat; 236ca3fa75bSLois Curfman McInnes PetscFunctionReturn(0); 237ca3fa75bSLois Curfman McInnes } 238ca3fa75bSLois Curfman McInnes 2398c778c55SBarry Smith PetscErrorCode MatDenseRestoreArray_MPIDense(Mat A,PetscScalar *array[]) 240ff14e315SSatish Balay { 24173a71a0fSBarry Smith Mat_MPIDense *a = (Mat_MPIDense*)A->data; 24273a71a0fSBarry Smith PetscErrorCode ierr; 24373a71a0fSBarry Smith 2443a40ed3dSBarry Smith PetscFunctionBegin; 2458c778c55SBarry Smith ierr = MatDenseRestoreArray(a->A,array);CHKERRQ(ierr); 2463a40ed3dSBarry Smith PetscFunctionReturn(0); 247ff14e315SSatish Balay } 248ff14e315SSatish Balay 249dfbe8321SBarry Smith PetscErrorCode MatAssemblyBegin_MPIDense(Mat mat,MatAssemblyType mode) 2508965ea79SLois Curfman McInnes { 25139ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 252ce94432eSBarry Smith MPI_Comm comm; 253dfbe8321SBarry Smith PetscErrorCode ierr; 25413f74950SBarry Smith PetscInt nstash,reallocs; 2558965ea79SLois Curfman McInnes InsertMode addv; 2568965ea79SLois Curfman McInnes 2573a40ed3dSBarry Smith PetscFunctionBegin; 258ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)mat,&comm);CHKERRQ(ierr); 2598965ea79SLois Curfman McInnes /* make sure all processors are either in INSERTMODE or ADDMODE */ 260b2566f29SBarry Smith ierr = MPIU_Allreduce((PetscEnum*)&mat->insertmode,(PetscEnum*)&addv,1,MPIU_ENUM,MPI_BOR,comm);CHKERRQ(ierr); 261ce94432eSBarry Smith if (addv == (ADD_VALUES|INSERT_VALUES)) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Cannot mix adds/inserts on different procs"); 262e0fa3b82SLois Curfman McInnes mat->insertmode = addv; /* in case this processor had no cache */ 2638965ea79SLois Curfman McInnes 264d0f46423SBarry Smith ierr = MatStashScatterBegin_Private(mat,&mat->stash,mat->rmap->range);CHKERRQ(ierr); 2658798bf22SSatish Balay ierr = MatStashGetInfo_Private(&mat->stash,&nstash,&reallocs);CHKERRQ(ierr); 266ae15b995SBarry Smith ierr = PetscInfo2(mdn->A,"Stash has %D entries, uses %D mallocs.\n",nstash,reallocs);CHKERRQ(ierr); 2673a40ed3dSBarry Smith PetscFunctionReturn(0); 2688965ea79SLois Curfman McInnes } 2698965ea79SLois Curfman McInnes 270dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_MPIDense(Mat mat,MatAssemblyType mode) 2718965ea79SLois Curfman McInnes { 27239ddd567SLois Curfman McInnes Mat_MPIDense *mdn=(Mat_MPIDense*)mat->data; 2736849ba73SBarry Smith PetscErrorCode ierr; 27413f74950SBarry Smith PetscInt i,*row,*col,flg,j,rstart,ncols; 27513f74950SBarry Smith PetscMPIInt n; 27687828ca2SBarry Smith PetscScalar *val; 2778965ea79SLois Curfman McInnes 2783a40ed3dSBarry Smith PetscFunctionBegin; 2798965ea79SLois Curfman McInnes /* wait on receives */ 2807ef1d9bdSSatish Balay while (1) { 2818798bf22SSatish Balay ierr = MatStashScatterGetMesg_Private(&mat->stash,&n,&row,&col,&val,&flg);CHKERRQ(ierr); 2827ef1d9bdSSatish Balay if (!flg) break; 2838965ea79SLois Curfman McInnes 2847ef1d9bdSSatish Balay for (i=0; i<n;) { 2857ef1d9bdSSatish Balay /* Now identify the consecutive vals belonging to the same row */ 2862205254eSKarl Rupp for (j=i,rstart=row[j]; j<n; j++) { 2872205254eSKarl Rupp if (row[j] != rstart) break; 2882205254eSKarl Rupp } 2897ef1d9bdSSatish Balay if (j < n) ncols = j-i; 2907ef1d9bdSSatish Balay else ncols = n-i; 2917ef1d9bdSSatish Balay /* Now assemble all these values with a single function call */ 2924b4eb8d3SJed Brown ierr = MatSetValues_MPIDense(mat,1,row+i,ncols,col+i,val+i,mat->insertmode);CHKERRQ(ierr); 2937ef1d9bdSSatish Balay i = j; 2948965ea79SLois Curfman McInnes } 2957ef1d9bdSSatish Balay } 2968798bf22SSatish Balay ierr = MatStashScatterEnd_Private(&mat->stash);CHKERRQ(ierr); 2978965ea79SLois Curfman McInnes 29839ddd567SLois Curfman McInnes ierr = MatAssemblyBegin(mdn->A,mode);CHKERRQ(ierr); 29939ddd567SLois Curfman McInnes ierr = MatAssemblyEnd(mdn->A,mode);CHKERRQ(ierr); 3008965ea79SLois Curfman McInnes 3016d4a8577SBarry Smith if (!mat->was_assembled && mode == MAT_FINAL_ASSEMBLY) { 30239ddd567SLois Curfman McInnes ierr = MatSetUpMultiply_MPIDense(mat);CHKERRQ(ierr); 3038965ea79SLois Curfman McInnes } 3043a40ed3dSBarry Smith PetscFunctionReturn(0); 3058965ea79SLois Curfman McInnes } 3068965ea79SLois Curfman McInnes 307dfbe8321SBarry Smith PetscErrorCode MatZeroEntries_MPIDense(Mat A) 3088965ea79SLois Curfman McInnes { 309dfbe8321SBarry Smith PetscErrorCode ierr; 31039ddd567SLois Curfman McInnes Mat_MPIDense *l = (Mat_MPIDense*)A->data; 3113a40ed3dSBarry Smith 3123a40ed3dSBarry Smith PetscFunctionBegin; 3133a40ed3dSBarry Smith ierr = MatZeroEntries(l->A);CHKERRQ(ierr); 3143a40ed3dSBarry Smith PetscFunctionReturn(0); 3158965ea79SLois Curfman McInnes } 3168965ea79SLois Curfman McInnes 3178965ea79SLois Curfman McInnes /* the code does not do the diagonal entries correctly unless the 3188965ea79SLois Curfman McInnes matrix is square and the column and row owerships are identical. 3198965ea79SLois Curfman McInnes This is a BUG. The only way to fix it seems to be to access 3203501a2bdSLois Curfman McInnes mdn->A and mdn->B directly and not through the MatZeroRows() 3218965ea79SLois Curfman McInnes routine. 3228965ea79SLois Curfman McInnes */ 3232b40b63fSBarry Smith PetscErrorCode MatZeroRows_MPIDense(Mat A,PetscInt N,const PetscInt rows[],PetscScalar diag,Vec x,Vec b) 3248965ea79SLois Curfman McInnes { 32539ddd567SLois Curfman McInnes Mat_MPIDense *l = (Mat_MPIDense*)A->data; 3266849ba73SBarry Smith PetscErrorCode ierr; 327d0f46423SBarry Smith PetscInt i,*owners = A->rmap->range; 32876ec1555SBarry Smith PetscInt *sizes,j,idx,nsends; 32913f74950SBarry Smith PetscInt nmax,*svalues,*starts,*owner,nrecvs; 3307adad957SLisandro Dalcin PetscInt *rvalues,tag = ((PetscObject)A)->tag,count,base,slen,*source; 33113f74950SBarry Smith PetscInt *lens,*lrows,*values; 33213f74950SBarry Smith PetscMPIInt n,imdex,rank = l->rank,size = l->size; 333ce94432eSBarry Smith MPI_Comm comm; 3348965ea79SLois Curfman McInnes MPI_Request *send_waits,*recv_waits; 3358965ea79SLois Curfman McInnes MPI_Status recv_status,*send_status; 336ace3abfcSBarry Smith PetscBool found; 33797b48c8fSBarry Smith const PetscScalar *xx; 33897b48c8fSBarry Smith PetscScalar *bb; 3398965ea79SLois Curfman McInnes 3403a40ed3dSBarry Smith PetscFunctionBegin; 341ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 342ce94432eSBarry Smith if (A->rmap->N != A->cmap->N) SETERRQ(comm,PETSC_ERR_SUP,"Only handles square matrices"); 343b9679d65SBarry Smith if (A->rmap->n != A->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only handles matrices with identical column and row ownership"); 3448965ea79SLois Curfman McInnes /* first count number of contributors to each processor */ 345f628708eSJed Brown ierr = PetscCalloc1(2*size,&sizes);CHKERRQ(ierr); 346037dbc42SBarry Smith ierr = PetscMalloc1(N+1,&owner);CHKERRQ(ierr); /* see note*/ 3478965ea79SLois Curfman McInnes for (i=0; i<N; i++) { 3488965ea79SLois Curfman McInnes idx = rows[i]; 34935d8aa7fSBarry Smith found = PETSC_FALSE; 3508965ea79SLois Curfman McInnes for (j=0; j<size; j++) { 3518965ea79SLois Curfman McInnes if (idx >= owners[j] && idx < owners[j+1]) { 35276ec1555SBarry Smith sizes[2*j]++; sizes[2*j+1] = 1; owner[i] = j; found = PETSC_TRUE; break; 3538965ea79SLois Curfman McInnes } 3548965ea79SLois Curfman McInnes } 355e32f2f54SBarry Smith if (!found) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Index out of range"); 3568965ea79SLois Curfman McInnes } 3572205254eSKarl Rupp nsends = 0; 35876ec1555SBarry Smith for (i=0; i<size; i++) nsends += sizes[2*i+1]; 3598965ea79SLois Curfman McInnes 3608965ea79SLois Curfman McInnes /* inform other processors of number of messages and max length*/ 36176ec1555SBarry Smith ierr = PetscMaxSum(comm,sizes,&nmax,&nrecvs);CHKERRQ(ierr); 3628965ea79SLois Curfman McInnes 3638965ea79SLois Curfman McInnes /* post receives: */ 364785e854fSJed Brown ierr = PetscMalloc1((nrecvs+1)*(nmax+1),&rvalues);CHKERRQ(ierr); 365854ce69bSBarry Smith ierr = PetscMalloc1(nrecvs+1,&recv_waits);CHKERRQ(ierr); 3668965ea79SLois Curfman McInnes for (i=0; i<nrecvs; i++) { 36713f74950SBarry Smith ierr = MPI_Irecv(rvalues+nmax*i,nmax,MPIU_INT,MPI_ANY_SOURCE,tag,comm,recv_waits+i);CHKERRQ(ierr); 3688965ea79SLois Curfman McInnes } 3698965ea79SLois Curfman McInnes 3708965ea79SLois Curfman McInnes /* do sends: 3718965ea79SLois Curfman McInnes 1) starts[i] gives the starting index in svalues for stuff going to 3728965ea79SLois Curfman McInnes the ith processor 3738965ea79SLois Curfman McInnes */ 374854ce69bSBarry Smith ierr = PetscMalloc1(N+1,&svalues);CHKERRQ(ierr); 375854ce69bSBarry Smith ierr = PetscMalloc1(nsends+1,&send_waits);CHKERRQ(ierr); 376854ce69bSBarry Smith ierr = PetscMalloc1(size+1,&starts);CHKERRQ(ierr); 3778965ea79SLois Curfman McInnes 3788965ea79SLois Curfman McInnes starts[0] = 0; 37976ec1555SBarry Smith for (i=1; i<size; i++) starts[i] = starts[i-1] + sizes[2*i-2]; 3802205254eSKarl Rupp for (i=0; i<N; i++) svalues[starts[owner[i]]++] = rows[i]; 3812205254eSKarl Rupp 3822205254eSKarl Rupp starts[0] = 0; 38376ec1555SBarry Smith for (i=1; i<size+1; i++) starts[i] = starts[i-1] + sizes[2*i-2]; 3848965ea79SLois Curfman McInnes count = 0; 3858965ea79SLois Curfman McInnes for (i=0; i<size; i++) { 38676ec1555SBarry Smith if (sizes[2*i+1]) { 38776ec1555SBarry Smith ierr = MPI_Isend(svalues+starts[i],sizes[2*i],MPIU_INT,i,tag,comm,send_waits+count++);CHKERRQ(ierr); 3888965ea79SLois Curfman McInnes } 3898965ea79SLois Curfman McInnes } 390606d414cSSatish Balay ierr = PetscFree(starts);CHKERRQ(ierr); 3918965ea79SLois Curfman McInnes 3928965ea79SLois Curfman McInnes base = owners[rank]; 3938965ea79SLois Curfman McInnes 3948965ea79SLois Curfman McInnes /* wait on receives */ 395dcca6d9dSJed Brown ierr = PetscMalloc2(nrecvs,&lens,nrecvs,&source);CHKERRQ(ierr); 39674ed9c26SBarry Smith count = nrecvs; 39774ed9c26SBarry Smith slen = 0; 3988965ea79SLois Curfman McInnes while (count) { 399ca161407SBarry Smith ierr = MPI_Waitany(nrecvs,recv_waits,&imdex,&recv_status);CHKERRQ(ierr); 4008965ea79SLois Curfman McInnes /* unpack receives into our local space */ 40113f74950SBarry Smith ierr = MPI_Get_count(&recv_status,MPIU_INT,&n);CHKERRQ(ierr); 4022205254eSKarl Rupp 4038965ea79SLois Curfman McInnes source[imdex] = recv_status.MPI_SOURCE; 4048965ea79SLois Curfman McInnes lens[imdex] = n; 4058965ea79SLois Curfman McInnes slen += n; 4068965ea79SLois Curfman McInnes count--; 4078965ea79SLois Curfman McInnes } 408606d414cSSatish Balay ierr = PetscFree(recv_waits);CHKERRQ(ierr); 4098965ea79SLois Curfman McInnes 4108965ea79SLois Curfman McInnes /* move the data into the send scatter */ 411854ce69bSBarry Smith ierr = PetscMalloc1(slen+1,&lrows);CHKERRQ(ierr); 4128965ea79SLois Curfman McInnes count = 0; 4138965ea79SLois Curfman McInnes for (i=0; i<nrecvs; i++) { 4148965ea79SLois Curfman McInnes values = rvalues + i*nmax; 4158965ea79SLois Curfman McInnes for (j=0; j<lens[i]; j++) { 4168965ea79SLois Curfman McInnes lrows[count++] = values[j] - base; 4178965ea79SLois Curfman McInnes } 4188965ea79SLois Curfman McInnes } 419606d414cSSatish Balay ierr = PetscFree(rvalues);CHKERRQ(ierr); 42074ed9c26SBarry Smith ierr = PetscFree2(lens,source);CHKERRQ(ierr); 421606d414cSSatish Balay ierr = PetscFree(owner);CHKERRQ(ierr); 42276ec1555SBarry Smith ierr = PetscFree(sizes);CHKERRQ(ierr); 4238965ea79SLois Curfman McInnes 42497b48c8fSBarry Smith /* fix right hand side if needed */ 42597b48c8fSBarry Smith if (x && b) { 42697b48c8fSBarry Smith ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr); 42797b48c8fSBarry Smith ierr = VecGetArray(b,&bb);CHKERRQ(ierr); 42897b48c8fSBarry Smith for (i=0; i<slen; i++) { 42997b48c8fSBarry Smith bb[lrows[i]] = diag*xx[lrows[i]]; 43097b48c8fSBarry Smith } 43197b48c8fSBarry Smith ierr = VecRestoreArrayRead(x,&xx);CHKERRQ(ierr); 43297b48c8fSBarry Smith ierr = VecRestoreArray(b,&bb);CHKERRQ(ierr); 43397b48c8fSBarry Smith } 43497b48c8fSBarry Smith 4358965ea79SLois Curfman McInnes /* actually zap the local rows */ 436b9679d65SBarry Smith ierr = MatZeroRows(l->A,slen,lrows,0.0,0,0);CHKERRQ(ierr); 437e2eb51b1SBarry Smith if (diag != 0.0) { 438b9679d65SBarry Smith Mat_SeqDense *ll = (Mat_SeqDense*)l->A->data; 439b9679d65SBarry Smith PetscInt m = ll->lda, i; 440b9679d65SBarry Smith 441b9679d65SBarry Smith for (i=0; i<slen; i++) { 442b9679d65SBarry Smith ll->v[lrows[i] + m*(A->cmap->rstart + lrows[i])] = diag; 443b9679d65SBarry Smith } 444b9679d65SBarry Smith } 445606d414cSSatish Balay ierr = PetscFree(lrows);CHKERRQ(ierr); 4468965ea79SLois Curfman McInnes 4478965ea79SLois Curfman McInnes /* wait on sends */ 4488965ea79SLois Curfman McInnes if (nsends) { 449785e854fSJed Brown ierr = PetscMalloc1(nsends,&send_status);CHKERRQ(ierr); 450ca161407SBarry Smith ierr = MPI_Waitall(nsends,send_waits,send_status);CHKERRQ(ierr); 451606d414cSSatish Balay ierr = PetscFree(send_status);CHKERRQ(ierr); 4528965ea79SLois Curfman McInnes } 453606d414cSSatish Balay ierr = PetscFree(send_waits);CHKERRQ(ierr); 454606d414cSSatish Balay ierr = PetscFree(svalues);CHKERRQ(ierr); 4553a40ed3dSBarry Smith PetscFunctionReturn(0); 4568965ea79SLois Curfman McInnes } 4578965ea79SLois Curfman McInnes 458cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatMult_SeqDense(Mat,Vec,Vec); 459cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatMultAdd_SeqDense(Mat,Vec,Vec,Vec); 460cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatMultTranspose_SeqDense(Mat,Vec,Vec); 461cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatMultTransposeAdd_SeqDense(Mat,Vec,Vec,Vec); 462cc2e6a90SBarry Smith 463dfbe8321SBarry Smith PetscErrorCode MatMult_MPIDense(Mat mat,Vec xx,Vec yy) 4648965ea79SLois Curfman McInnes { 46539ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 466dfbe8321SBarry Smith PetscErrorCode ierr; 467c456f294SBarry Smith 4683a40ed3dSBarry Smith PetscFunctionBegin; 469ca9f406cSSatish Balay ierr = VecScatterBegin(mdn->Mvctx,xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 470ca9f406cSSatish Balay ierr = VecScatterEnd(mdn->Mvctx,xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 471cc2e6a90SBarry Smith ierr = MatMult_SeqDense(mdn->A,mdn->lvec,yy);CHKERRQ(ierr); 4723a40ed3dSBarry Smith PetscFunctionReturn(0); 4738965ea79SLois Curfman McInnes } 4748965ea79SLois Curfman McInnes 475dfbe8321SBarry Smith PetscErrorCode MatMultAdd_MPIDense(Mat mat,Vec xx,Vec yy,Vec zz) 4768965ea79SLois Curfman McInnes { 47739ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 478dfbe8321SBarry Smith PetscErrorCode ierr; 479c456f294SBarry Smith 4803a40ed3dSBarry Smith PetscFunctionBegin; 481ca9f406cSSatish Balay ierr = VecScatterBegin(mdn->Mvctx,xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 482ca9f406cSSatish Balay ierr = VecScatterEnd(mdn->Mvctx,xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 483cc2e6a90SBarry Smith ierr = MatMultAdd_SeqDense(mdn->A,mdn->lvec,yy,zz);CHKERRQ(ierr); 4843a40ed3dSBarry Smith PetscFunctionReturn(0); 4858965ea79SLois Curfman McInnes } 4868965ea79SLois Curfman McInnes 487dfbe8321SBarry Smith PetscErrorCode MatMultTranspose_MPIDense(Mat A,Vec xx,Vec yy) 488096963f5SLois Curfman McInnes { 489096963f5SLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense*)A->data; 490dfbe8321SBarry Smith PetscErrorCode ierr; 49187828ca2SBarry Smith PetscScalar zero = 0.0; 492096963f5SLois Curfman McInnes 4933a40ed3dSBarry Smith PetscFunctionBegin; 4942dcb1b2aSMatthew Knepley ierr = VecSet(yy,zero);CHKERRQ(ierr); 495cc2e6a90SBarry Smith ierr = MatMultTranspose_SeqDense(a->A,xx,a->lvec);CHKERRQ(ierr); 496ca9f406cSSatish Balay ierr = VecScatterBegin(a->Mvctx,a->lvec,yy,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 497ca9f406cSSatish Balay ierr = VecScatterEnd(a->Mvctx,a->lvec,yy,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 4983a40ed3dSBarry Smith PetscFunctionReturn(0); 499096963f5SLois Curfman McInnes } 500096963f5SLois Curfman McInnes 501dfbe8321SBarry Smith PetscErrorCode MatMultTransposeAdd_MPIDense(Mat A,Vec xx,Vec yy,Vec zz) 502096963f5SLois Curfman McInnes { 503096963f5SLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense*)A->data; 504dfbe8321SBarry Smith PetscErrorCode ierr; 505096963f5SLois Curfman McInnes 5063a40ed3dSBarry Smith PetscFunctionBegin; 5073501a2bdSLois Curfman McInnes ierr = VecCopy(yy,zz);CHKERRQ(ierr); 508cc2e6a90SBarry Smith ierr = MatMultTranspose_SeqDense(a->A,xx,a->lvec);CHKERRQ(ierr); 509ca9f406cSSatish Balay ierr = VecScatterBegin(a->Mvctx,a->lvec,zz,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 510ca9f406cSSatish Balay ierr = VecScatterEnd(a->Mvctx,a->lvec,zz,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 5113a40ed3dSBarry Smith PetscFunctionReturn(0); 512096963f5SLois Curfman McInnes } 513096963f5SLois Curfman McInnes 514dfbe8321SBarry Smith PetscErrorCode MatGetDiagonal_MPIDense(Mat A,Vec v) 5158965ea79SLois Curfman McInnes { 51639ddd567SLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense*)A->data; 517096963f5SLois Curfman McInnes Mat_SeqDense *aloc = (Mat_SeqDense*)a->A->data; 518dfbe8321SBarry Smith PetscErrorCode ierr; 519d0f46423SBarry Smith PetscInt len,i,n,m = A->rmap->n,radd; 52087828ca2SBarry Smith PetscScalar *x,zero = 0.0; 521ed3cc1f0SBarry Smith 5223a40ed3dSBarry Smith PetscFunctionBegin; 5232dcb1b2aSMatthew Knepley ierr = VecSet(v,zero);CHKERRQ(ierr); 5241ebc52fbSHong Zhang ierr = VecGetArray(v,&x);CHKERRQ(ierr); 525096963f5SLois Curfman McInnes ierr = VecGetSize(v,&n);CHKERRQ(ierr); 526e32f2f54SBarry Smith if (n != A->rmap->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec"); 527d0f46423SBarry Smith len = PetscMin(a->A->rmap->n,a->A->cmap->n); 528d0f46423SBarry Smith radd = A->rmap->rstart*m; 52944cd7ae7SLois Curfman McInnes for (i=0; i<len; i++) { 530096963f5SLois Curfman McInnes x[i] = aloc->v[radd + i*m + i]; 531096963f5SLois Curfman McInnes } 5321ebc52fbSHong Zhang ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 5333a40ed3dSBarry Smith PetscFunctionReturn(0); 5348965ea79SLois Curfman McInnes } 5358965ea79SLois Curfman McInnes 536dfbe8321SBarry Smith PetscErrorCode MatDestroy_MPIDense(Mat mat) 5378965ea79SLois Curfman McInnes { 5383501a2bdSLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 539dfbe8321SBarry Smith PetscErrorCode ierr; 540ed3cc1f0SBarry Smith 5413a40ed3dSBarry Smith PetscFunctionBegin; 542aa482453SBarry Smith #if defined(PETSC_USE_LOG) 543d0f46423SBarry Smith PetscLogObjectState((PetscObject)mat,"Rows=%D, Cols=%D",mat->rmap->N,mat->cmap->N); 5448965ea79SLois Curfman McInnes #endif 5458798bf22SSatish Balay ierr = MatStashDestroy_Private(&mat->stash);CHKERRQ(ierr); 5466bf464f9SBarry Smith ierr = MatDestroy(&mdn->A);CHKERRQ(ierr); 5476bf464f9SBarry Smith ierr = VecDestroy(&mdn->lvec);CHKERRQ(ierr); 5486bf464f9SBarry Smith ierr = VecScatterDestroy(&mdn->Mvctx);CHKERRQ(ierr); 54901b82886SBarry Smith 550bf0cc555SLisandro Dalcin ierr = PetscFree(mat->data);CHKERRQ(ierr); 551dbd8c25aSHong Zhang ierr = PetscObjectChangeTypeName((PetscObject)mat,0);CHKERRQ(ierr); 5528baccfbdSHong Zhang 5538baccfbdSHong Zhang ierr = PetscObjectComposeFunction((PetscObject)mat,"MatDenseGetArray_C",NULL);CHKERRQ(ierr); 554d3042a70SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatDensePlaceArray_C",NULL);CHKERRQ(ierr); 555d3042a70SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatDenseResetArray_C",NULL);CHKERRQ(ierr); 556d3042a70SBarry Smith 5578baccfbdSHong Zhang ierr = PetscObjectComposeFunction((PetscObject)mat,"MatDenseRestoreArray_C",NULL);CHKERRQ(ierr); 5588baccfbdSHong Zhang #if defined(PETSC_HAVE_ELEMENTAL) 5598baccfbdSHong Zhang ierr = PetscObjectComposeFunction((PetscObject)mat,"MatConvert_mpidense_elemental_C",NULL);CHKERRQ(ierr); 5608baccfbdSHong Zhang #endif 561bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMPIDenseSetPreallocation_C",NULL);CHKERRQ(ierr); 562bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMatMult_mpiaij_mpidense_C",NULL);CHKERRQ(ierr); 563bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMatMultSymbolic_mpiaij_mpidense_C",NULL);CHKERRQ(ierr); 564bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMatMultNumeric_mpiaij_mpidense_C",NULL);CHKERRQ(ierr); 5658baccfbdSHong Zhang ierr = PetscObjectComposeFunction((PetscObject)mat,"MatTransposeMatMult_mpiaij_mpidense_C",NULL);CHKERRQ(ierr); 5668baccfbdSHong Zhang ierr = PetscObjectComposeFunction((PetscObject)mat,"MatTransposeMatMultSymbolic_mpiaij_mpidense_C",NULL);CHKERRQ(ierr); 5678baccfbdSHong Zhang ierr = PetscObjectComposeFunction((PetscObject)mat,"MatTransposeMatMultNumeric_mpiaij_mpidense_C",NULL);CHKERRQ(ierr); 568*86aefd0dSHong Zhang ierr = PetscObjectComposeFunction((PetscObject)mat,"MatDenseGetColumn_C",NULL);CHKERRQ(ierr); 569*86aefd0dSHong Zhang ierr = PetscObjectComposeFunction((PetscObject)mat,"MatDenseRestoreColumn_C",NULL);CHKERRQ(ierr); 5703a40ed3dSBarry Smith PetscFunctionReturn(0); 5718965ea79SLois Curfman McInnes } 57239ddd567SLois Curfman McInnes 5736849ba73SBarry Smith static PetscErrorCode MatView_MPIDense_Binary(Mat mat,PetscViewer viewer) 5748965ea79SLois Curfman McInnes { 57539ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 576dfbe8321SBarry Smith PetscErrorCode ierr; 577aa05aa95SBarry Smith PetscViewerFormat format; 578aa05aa95SBarry Smith int fd; 579d0f46423SBarry Smith PetscInt header[4],mmax,N = mat->cmap->N,i,j,m,k; 580aa05aa95SBarry Smith PetscMPIInt rank,tag = ((PetscObject)viewer)->tag,size; 581578230a0SSatish Balay PetscScalar *work,*v,*vv; 582aa05aa95SBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)mdn->A->data; 5837056b6fcSBarry Smith 5843a40ed3dSBarry Smith PetscFunctionBegin; 58539ddd567SLois Curfman McInnes if (mdn->size == 1) { 58639ddd567SLois Curfman McInnes ierr = MatView(mdn->A,viewer);CHKERRQ(ierr); 587aa05aa95SBarry Smith } else { 588aa05aa95SBarry Smith ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 589ce94432eSBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)mat),&rank);CHKERRQ(ierr); 590ce94432eSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)mat),&size);CHKERRQ(ierr); 591aa05aa95SBarry Smith 592aa05aa95SBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 593f4403165SShri Abhyankar if (format == PETSC_VIEWER_NATIVE) { 594aa05aa95SBarry Smith 595aa05aa95SBarry Smith if (!rank) { 596aa05aa95SBarry Smith /* store the matrix as a dense matrix */ 5970700a824SBarry Smith header[0] = MAT_FILE_CLASSID; 598d0f46423SBarry Smith header[1] = mat->rmap->N; 599aa05aa95SBarry Smith header[2] = N; 600aa05aa95SBarry Smith header[3] = MATRIX_BINARY_FORMAT_DENSE; 601aa05aa95SBarry Smith ierr = PetscBinaryWrite(fd,header,4,PETSC_INT,PETSC_TRUE);CHKERRQ(ierr); 602aa05aa95SBarry Smith 603aa05aa95SBarry Smith /* get largest work array needed for transposing array */ 604d0f46423SBarry Smith mmax = mat->rmap->n; 605aa05aa95SBarry Smith for (i=1; i<size; i++) { 606d0f46423SBarry Smith mmax = PetscMax(mmax,mat->rmap->range[i+1] - mat->rmap->range[i]); 6078965ea79SLois Curfman McInnes } 608785e854fSJed Brown ierr = PetscMalloc1(mmax*N,&work);CHKERRQ(ierr); 609aa05aa95SBarry Smith 610aa05aa95SBarry Smith /* write out local array, by rows */ 611d0f46423SBarry Smith m = mat->rmap->n; 612aa05aa95SBarry Smith v = a->v; 613aa05aa95SBarry Smith for (j=0; j<N; j++) { 614aa05aa95SBarry Smith for (i=0; i<m; i++) { 615578230a0SSatish Balay work[j + i*N] = *v++; 616aa05aa95SBarry Smith } 617aa05aa95SBarry Smith } 618aa05aa95SBarry Smith ierr = PetscBinaryWrite(fd,work,m*N,PETSC_SCALAR,PETSC_FALSE);CHKERRQ(ierr); 619aa05aa95SBarry Smith /* get largest work array to receive messages from other processes, excludes process zero */ 620aa05aa95SBarry Smith mmax = 0; 621aa05aa95SBarry Smith for (i=1; i<size; i++) { 622d0f46423SBarry Smith mmax = PetscMax(mmax,mat->rmap->range[i+1] - mat->rmap->range[i]); 623aa05aa95SBarry Smith } 624785e854fSJed Brown ierr = PetscMalloc1(mmax*N,&vv);CHKERRQ(ierr); 625aa05aa95SBarry Smith for (k = 1; k < size; k++) { 626f8009846SMatthew Knepley v = vv; 627d0f46423SBarry Smith m = mat->rmap->range[k+1] - mat->rmap->range[k]; 628ce94432eSBarry Smith ierr = MPIULong_Recv(v,m*N,MPIU_SCALAR,k,tag,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr); 629aa05aa95SBarry Smith 630aa05aa95SBarry Smith for (j = 0; j < N; j++) { 631aa05aa95SBarry Smith for (i = 0; i < m; i++) { 632578230a0SSatish Balay work[j + i*N] = *v++; 633aa05aa95SBarry Smith } 634aa05aa95SBarry Smith } 635aa05aa95SBarry Smith ierr = PetscBinaryWrite(fd,work,m*N,PETSC_SCALAR,PETSC_FALSE);CHKERRQ(ierr); 636aa05aa95SBarry Smith } 637aa05aa95SBarry Smith ierr = PetscFree(work);CHKERRQ(ierr); 638578230a0SSatish Balay ierr = PetscFree(vv);CHKERRQ(ierr); 639aa05aa95SBarry Smith } else { 640ce94432eSBarry Smith ierr = MPIULong_Send(a->v,mat->rmap->n*mat->cmap->N,MPIU_SCALAR,0,tag,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr); 641aa05aa95SBarry Smith } 6426a9046bcSBarry Smith } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"To store a parallel dense matrix you must first call PetscViewerPushFormat(viewer,PETSC_VIEWER_NATIVE)"); 643aa05aa95SBarry Smith } 6443a40ed3dSBarry Smith PetscFunctionReturn(0); 6458965ea79SLois Curfman McInnes } 6468965ea79SLois Curfman McInnes 6477da1fb6eSBarry Smith extern PetscErrorCode MatView_SeqDense(Mat,PetscViewer); 6489804daf3SBarry Smith #include <petscdraw.h> 6496849ba73SBarry Smith static PetscErrorCode MatView_MPIDense_ASCIIorDraworSocket(Mat mat,PetscViewer viewer) 6508965ea79SLois Curfman McInnes { 65139ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 652dfbe8321SBarry Smith PetscErrorCode ierr; 6537da1fb6eSBarry Smith PetscMPIInt rank = mdn->rank; 65419fd82e9SBarry Smith PetscViewerType vtype; 655ace3abfcSBarry Smith PetscBool iascii,isdraw; 656b0a32e0cSBarry Smith PetscViewer sviewer; 657f3ef73ceSBarry Smith PetscViewerFormat format; 6588965ea79SLois Curfman McInnes 6593a40ed3dSBarry Smith PetscFunctionBegin; 660251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 661251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr); 66232077d6dSBarry Smith if (iascii) { 663b0a32e0cSBarry Smith ierr = PetscViewerGetType(viewer,&vtype);CHKERRQ(ierr); 664b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 665456192e2SBarry Smith if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) { 6664e220ebcSLois Curfman McInnes MatInfo info; 667888f2ed8SSatish Balay ierr = MatGetInfo(mat,MAT_LOCAL,&info);CHKERRQ(ierr); 6681575c14dSBarry Smith ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr); 6697b23a99aSBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer," [%d] local rows %D nz %D nz alloced %D mem %D \n",rank,mat->rmap->n,(PetscInt)info.nz_used,(PetscInt)info.nz_allocated,(PetscInt)info.memory);CHKERRQ(ierr); 670b0a32e0cSBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 6711575c14dSBarry Smith ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr); 6725d00a290SHong Zhang ierr = VecScatterView(mdn->Mvctx,viewer);CHKERRQ(ierr); 6733a40ed3dSBarry Smith PetscFunctionReturn(0); 674fb9695e5SSatish Balay } else if (format == PETSC_VIEWER_ASCII_INFO) { 6753a40ed3dSBarry Smith PetscFunctionReturn(0); 6768965ea79SLois Curfman McInnes } 677f1af5d2fSBarry Smith } else if (isdraw) { 678b0a32e0cSBarry Smith PetscDraw draw; 679ace3abfcSBarry Smith PetscBool isnull; 680f1af5d2fSBarry Smith 681b0a32e0cSBarry Smith ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 682b0a32e0cSBarry Smith ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); 683f1af5d2fSBarry Smith if (isnull) PetscFunctionReturn(0); 684f1af5d2fSBarry Smith } 68577ed5343SBarry Smith 6867da1fb6eSBarry Smith { 6878965ea79SLois Curfman McInnes /* assemble the entire matrix onto first processor. */ 6888965ea79SLois Curfman McInnes Mat A; 689d0f46423SBarry Smith PetscInt M = mat->rmap->N,N = mat->cmap->N,m,row,i,nz; 690ba8c8a56SBarry Smith PetscInt *cols; 691ba8c8a56SBarry Smith PetscScalar *vals; 6928965ea79SLois Curfman McInnes 693ce94432eSBarry Smith ierr = MatCreate(PetscObjectComm((PetscObject)mat),&A);CHKERRQ(ierr); 6948965ea79SLois Curfman McInnes if (!rank) { 695f69a0ea3SMatthew Knepley ierr = MatSetSizes(A,M,N,M,N);CHKERRQ(ierr); 6963a40ed3dSBarry Smith } else { 697f69a0ea3SMatthew Knepley ierr = MatSetSizes(A,0,0,M,N);CHKERRQ(ierr); 6988965ea79SLois Curfman McInnes } 6997adad957SLisandro Dalcin /* Since this is a temporary matrix, MATMPIDENSE instead of ((PetscObject)A)->type_name here is probably acceptable. */ 700878740d9SKris Buschelman ierr = MatSetType(A,MATMPIDENSE);CHKERRQ(ierr); 7010298fd71SBarry Smith ierr = MatMPIDenseSetPreallocation(A,NULL);CHKERRQ(ierr); 7023bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)A);CHKERRQ(ierr); 7038965ea79SLois Curfman McInnes 70439ddd567SLois Curfman McInnes /* Copy the matrix ... This isn't the most efficient means, 70539ddd567SLois Curfman McInnes but it's quick for now */ 70651022da4SBarry Smith A->insertmode = INSERT_VALUES; 7072205254eSKarl Rupp 7082205254eSKarl Rupp row = mat->rmap->rstart; 7092205254eSKarl Rupp m = mdn->A->rmap->n; 7108965ea79SLois Curfman McInnes for (i=0; i<m; i++) { 711ba8c8a56SBarry Smith ierr = MatGetRow_MPIDense(mat,row,&nz,&cols,&vals);CHKERRQ(ierr); 712ba8c8a56SBarry Smith ierr = MatSetValues_MPIDense(A,1,&row,nz,cols,vals,INSERT_VALUES);CHKERRQ(ierr); 713ba8c8a56SBarry Smith ierr = MatRestoreRow_MPIDense(mat,row,&nz,&cols,&vals);CHKERRQ(ierr); 71439ddd567SLois Curfman McInnes row++; 7158965ea79SLois Curfman McInnes } 7168965ea79SLois Curfman McInnes 7176d4a8577SBarry Smith ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 7186d4a8577SBarry Smith ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 7193f08860eSBarry Smith ierr = PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr); 720b9b97703SBarry Smith if (!rank) { 7211a9d3c3cSBarry Smith ierr = PetscObjectSetName((PetscObject)((Mat_MPIDense*)(A->data))->A,((PetscObject)mat)->name);CHKERRQ(ierr); 7227da1fb6eSBarry Smith ierr = MatView_SeqDense(((Mat_MPIDense*)(A->data))->A,sviewer);CHKERRQ(ierr); 7238965ea79SLois Curfman McInnes } 7243f08860eSBarry Smith ierr = PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr); 725b0a32e0cSBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 7266bf464f9SBarry Smith ierr = MatDestroy(&A);CHKERRQ(ierr); 7278965ea79SLois Curfman McInnes } 7283a40ed3dSBarry Smith PetscFunctionReturn(0); 7298965ea79SLois Curfman McInnes } 7308965ea79SLois Curfman McInnes 731dfbe8321SBarry Smith PetscErrorCode MatView_MPIDense(Mat mat,PetscViewer viewer) 7328965ea79SLois Curfman McInnes { 733dfbe8321SBarry Smith PetscErrorCode ierr; 734ace3abfcSBarry Smith PetscBool iascii,isbinary,isdraw,issocket; 7358965ea79SLois Curfman McInnes 736433994e6SBarry Smith PetscFunctionBegin; 737251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 738251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr); 739251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSOCKET,&issocket);CHKERRQ(ierr); 740251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr); 7410f5bd95cSBarry Smith 74232077d6dSBarry Smith if (iascii || issocket || isdraw) { 743f1af5d2fSBarry Smith ierr = MatView_MPIDense_ASCIIorDraworSocket(mat,viewer);CHKERRQ(ierr); 7440f5bd95cSBarry Smith } else if (isbinary) { 7453a40ed3dSBarry Smith ierr = MatView_MPIDense_Binary(mat,viewer);CHKERRQ(ierr); 74611aeaf0aSBarry Smith } 7473a40ed3dSBarry Smith PetscFunctionReturn(0); 7488965ea79SLois Curfman McInnes } 7498965ea79SLois Curfman McInnes 750dfbe8321SBarry Smith PetscErrorCode MatGetInfo_MPIDense(Mat A,MatInfoType flag,MatInfo *info) 7518965ea79SLois Curfman McInnes { 7523501a2bdSLois Curfman McInnes Mat_MPIDense *mat = (Mat_MPIDense*)A->data; 7533501a2bdSLois Curfman McInnes Mat mdn = mat->A; 754dfbe8321SBarry Smith PetscErrorCode ierr; 755329f5518SBarry Smith PetscReal isend[5],irecv[5]; 7568965ea79SLois Curfman McInnes 7573a40ed3dSBarry Smith PetscFunctionBegin; 7584e220ebcSLois Curfman McInnes info->block_size = 1.0; 7592205254eSKarl Rupp 7604e220ebcSLois Curfman McInnes ierr = MatGetInfo(mdn,MAT_LOCAL,info);CHKERRQ(ierr); 7612205254eSKarl Rupp 7624e220ebcSLois Curfman McInnes isend[0] = info->nz_used; isend[1] = info->nz_allocated; isend[2] = info->nz_unneeded; 7634e220ebcSLois Curfman McInnes isend[3] = info->memory; isend[4] = info->mallocs; 7648965ea79SLois Curfman McInnes if (flag == MAT_LOCAL) { 7654e220ebcSLois Curfman McInnes info->nz_used = isend[0]; 7664e220ebcSLois Curfman McInnes info->nz_allocated = isend[1]; 7674e220ebcSLois Curfman McInnes info->nz_unneeded = isend[2]; 7684e220ebcSLois Curfman McInnes info->memory = isend[3]; 7694e220ebcSLois Curfman McInnes info->mallocs = isend[4]; 7708965ea79SLois Curfman McInnes } else if (flag == MAT_GLOBAL_MAX) { 771b2566f29SBarry Smith ierr = MPIU_Allreduce(isend,irecv,5,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 7722205254eSKarl Rupp 7734e220ebcSLois Curfman McInnes info->nz_used = irecv[0]; 7744e220ebcSLois Curfman McInnes info->nz_allocated = irecv[1]; 7754e220ebcSLois Curfman McInnes info->nz_unneeded = irecv[2]; 7764e220ebcSLois Curfman McInnes info->memory = irecv[3]; 7774e220ebcSLois Curfman McInnes info->mallocs = irecv[4]; 7788965ea79SLois Curfman McInnes } else if (flag == MAT_GLOBAL_SUM) { 779b2566f29SBarry Smith ierr = MPIU_Allreduce(isend,irecv,5,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 7802205254eSKarl Rupp 7814e220ebcSLois Curfman McInnes info->nz_used = irecv[0]; 7824e220ebcSLois Curfman McInnes info->nz_allocated = irecv[1]; 7834e220ebcSLois Curfman McInnes info->nz_unneeded = irecv[2]; 7844e220ebcSLois Curfman McInnes info->memory = irecv[3]; 7854e220ebcSLois Curfman McInnes info->mallocs = irecv[4]; 7868965ea79SLois Curfman McInnes } 7874e220ebcSLois Curfman McInnes info->fill_ratio_given = 0; /* no parallel LU/ILU/Cholesky */ 7884e220ebcSLois Curfman McInnes info->fill_ratio_needed = 0; 7894e220ebcSLois Curfman McInnes info->factor_mallocs = 0; 7903a40ed3dSBarry Smith PetscFunctionReturn(0); 7918965ea79SLois Curfman McInnes } 7928965ea79SLois Curfman McInnes 793ace3abfcSBarry Smith PetscErrorCode MatSetOption_MPIDense(Mat A,MatOption op,PetscBool flg) 7948965ea79SLois Curfman McInnes { 79539ddd567SLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense*)A->data; 796dfbe8321SBarry Smith PetscErrorCode ierr; 7978965ea79SLois Curfman McInnes 7983a40ed3dSBarry Smith PetscFunctionBegin; 79912c028f9SKris Buschelman switch (op) { 800512a5fc5SBarry Smith case MAT_NEW_NONZERO_LOCATIONS: 80112c028f9SKris Buschelman case MAT_NEW_NONZERO_LOCATION_ERR: 80212c028f9SKris Buschelman case MAT_NEW_NONZERO_ALLOCATION_ERR: 80343674050SBarry Smith MatCheckPreallocated(A,1); 8044e0d8c25SBarry Smith ierr = MatSetOption(a->A,op,flg);CHKERRQ(ierr); 80512c028f9SKris Buschelman break; 80612c028f9SKris Buschelman case MAT_ROW_ORIENTED: 80743674050SBarry Smith MatCheckPreallocated(A,1); 8084e0d8c25SBarry Smith a->roworiented = flg; 8094e0d8c25SBarry Smith ierr = MatSetOption(a->A,op,flg);CHKERRQ(ierr); 81012c028f9SKris Buschelman break; 8114e0d8c25SBarry Smith case MAT_NEW_DIAGONALS: 81213fa8e87SLisandro Dalcin case MAT_KEEP_NONZERO_PATTERN: 81312c028f9SKris Buschelman case MAT_USE_HASH_TABLE: 814290bbb0aSBarry Smith ierr = PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);CHKERRQ(ierr); 81512c028f9SKris Buschelman break; 81612c028f9SKris Buschelman case MAT_IGNORE_OFF_PROC_ENTRIES: 8174e0d8c25SBarry Smith a->donotstash = flg; 81812c028f9SKris Buschelman break; 81977e54ba9SKris Buschelman case MAT_SYMMETRIC: 82077e54ba9SKris Buschelman case MAT_STRUCTURALLY_SYMMETRIC: 8219a4540c5SBarry Smith case MAT_HERMITIAN: 8229a4540c5SBarry Smith case MAT_SYMMETRY_ETERNAL: 823600fe468SBarry Smith case MAT_IGNORE_LOWER_TRIANGULAR: 824290bbb0aSBarry Smith ierr = PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);CHKERRQ(ierr); 82577e54ba9SKris Buschelman break; 82612c028f9SKris Buschelman default: 827e32f2f54SBarry Smith SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"unknown option %s",MatOptions[op]); 8283a40ed3dSBarry Smith } 8293a40ed3dSBarry Smith PetscFunctionReturn(0); 8308965ea79SLois Curfman McInnes } 8318965ea79SLois Curfman McInnes 8328965ea79SLois Curfman McInnes 833dfbe8321SBarry Smith PetscErrorCode MatDiagonalScale_MPIDense(Mat A,Vec ll,Vec rr) 8345b2fa520SLois Curfman McInnes { 8355b2fa520SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)A->data; 8365b2fa520SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)mdn->A->data; 837bca11509SBarry Smith const PetscScalar *l,*r; 838bca11509SBarry Smith PetscScalar x,*v; 839dfbe8321SBarry Smith PetscErrorCode ierr; 840d0f46423SBarry Smith PetscInt i,j,s2a,s3a,s2,s3,m=mdn->A->rmap->n,n=mdn->A->cmap->n; 8415b2fa520SLois Curfman McInnes 8425b2fa520SLois Curfman McInnes PetscFunctionBegin; 84372d926a5SLois Curfman McInnes ierr = MatGetLocalSize(A,&s2,&s3);CHKERRQ(ierr); 8445b2fa520SLois Curfman McInnes if (ll) { 84572d926a5SLois Curfman McInnes ierr = VecGetLocalSize(ll,&s2a);CHKERRQ(ierr); 846e32f2f54SBarry Smith if (s2a != s2) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Left scaling vector non-conforming local size, %d != %d.", s2a, s2); 847bca11509SBarry Smith ierr = VecGetArrayRead(ll,&l);CHKERRQ(ierr); 8485b2fa520SLois Curfman McInnes for (i=0; i<m; i++) { 8495b2fa520SLois Curfman McInnes x = l[i]; 8505b2fa520SLois Curfman McInnes v = mat->v + i; 8515b2fa520SLois Curfman McInnes for (j=0; j<n; j++) { (*v) *= x; v+= m;} 8525b2fa520SLois Curfman McInnes } 853bca11509SBarry Smith ierr = VecRestoreArrayRead(ll,&l);CHKERRQ(ierr); 854efee365bSSatish Balay ierr = PetscLogFlops(n*m);CHKERRQ(ierr); 8555b2fa520SLois Curfman McInnes } 8565b2fa520SLois Curfman McInnes if (rr) { 857175be7b4SMatthew Knepley ierr = VecGetLocalSize(rr,&s3a);CHKERRQ(ierr); 858e32f2f54SBarry Smith if (s3a != s3) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Right scaling vec non-conforming local size, %d != %d.", s3a, s3); 859ca9f406cSSatish Balay ierr = VecScatterBegin(mdn->Mvctx,rr,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 860ca9f406cSSatish Balay ierr = VecScatterEnd(mdn->Mvctx,rr,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 861bca11509SBarry Smith ierr = VecGetArrayRead(mdn->lvec,&r);CHKERRQ(ierr); 8625b2fa520SLois Curfman McInnes for (i=0; i<n; i++) { 8635b2fa520SLois Curfman McInnes x = r[i]; 8645b2fa520SLois Curfman McInnes v = mat->v + i*m; 8652205254eSKarl Rupp for (j=0; j<m; j++) (*v++) *= x; 8665b2fa520SLois Curfman McInnes } 867bca11509SBarry Smith ierr = VecRestoreArrayRead(mdn->lvec,&r);CHKERRQ(ierr); 868efee365bSSatish Balay ierr = PetscLogFlops(n*m);CHKERRQ(ierr); 8695b2fa520SLois Curfman McInnes } 8705b2fa520SLois Curfman McInnes PetscFunctionReturn(0); 8715b2fa520SLois Curfman McInnes } 8725b2fa520SLois Curfman McInnes 873dfbe8321SBarry Smith PetscErrorCode MatNorm_MPIDense(Mat A,NormType type,PetscReal *nrm) 874096963f5SLois Curfman McInnes { 8753501a2bdSLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)A->data; 8763501a2bdSLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)mdn->A->data; 877dfbe8321SBarry Smith PetscErrorCode ierr; 87813f74950SBarry Smith PetscInt i,j; 879329f5518SBarry Smith PetscReal sum = 0.0; 88087828ca2SBarry Smith PetscScalar *v = mat->v; 8813501a2bdSLois Curfman McInnes 8823a40ed3dSBarry Smith PetscFunctionBegin; 8833501a2bdSLois Curfman McInnes if (mdn->size == 1) { 884064f8208SBarry Smith ierr = MatNorm(mdn->A,type,nrm);CHKERRQ(ierr); 8853501a2bdSLois Curfman McInnes } else { 8863501a2bdSLois Curfman McInnes if (type == NORM_FROBENIUS) { 887d0f46423SBarry Smith for (i=0; i<mdn->A->cmap->n*mdn->A->rmap->n; i++) { 888329f5518SBarry Smith sum += PetscRealPart(PetscConj(*v)*(*v)); v++; 8893501a2bdSLois Curfman McInnes } 890b2566f29SBarry Smith ierr = MPIU_Allreduce(&sum,nrm,1,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 8918f1a2a5eSBarry Smith *nrm = PetscSqrtReal(*nrm); 892dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*mdn->A->cmap->n*mdn->A->rmap->n);CHKERRQ(ierr); 8933a40ed3dSBarry Smith } else if (type == NORM_1) { 894329f5518SBarry Smith PetscReal *tmp,*tmp2; 895dcca6d9dSJed Brown ierr = PetscMalloc2(A->cmap->N,&tmp,A->cmap->N,&tmp2);CHKERRQ(ierr); 89674ed9c26SBarry Smith ierr = PetscMemzero(tmp,A->cmap->N*sizeof(PetscReal));CHKERRQ(ierr); 89774ed9c26SBarry Smith ierr = PetscMemzero(tmp2,A->cmap->N*sizeof(PetscReal));CHKERRQ(ierr); 898064f8208SBarry Smith *nrm = 0.0; 8993501a2bdSLois Curfman McInnes v = mat->v; 900d0f46423SBarry Smith for (j=0; j<mdn->A->cmap->n; j++) { 901d0f46423SBarry Smith for (i=0; i<mdn->A->rmap->n; i++) { 90267e560aaSBarry Smith tmp[j] += PetscAbsScalar(*v); v++; 9033501a2bdSLois Curfman McInnes } 9043501a2bdSLois Curfman McInnes } 905b2566f29SBarry Smith ierr = MPIU_Allreduce(tmp,tmp2,A->cmap->N,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 906d0f46423SBarry Smith for (j=0; j<A->cmap->N; j++) { 907064f8208SBarry Smith if (tmp2[j] > *nrm) *nrm = tmp2[j]; 9083501a2bdSLois Curfman McInnes } 9098627564fSBarry Smith ierr = PetscFree2(tmp,tmp2);CHKERRQ(ierr); 910d0f46423SBarry Smith ierr = PetscLogFlops(A->cmap->n*A->rmap->n);CHKERRQ(ierr); 9113a40ed3dSBarry Smith } else if (type == NORM_INFINITY) { /* max row norm */ 912329f5518SBarry Smith PetscReal ntemp; 9133501a2bdSLois Curfman McInnes ierr = MatNorm(mdn->A,type,&ntemp);CHKERRQ(ierr); 914b2566f29SBarry Smith ierr = MPIU_Allreduce(&ntemp,nrm,1,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 915ce94432eSBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support for two norm"); 9163501a2bdSLois Curfman McInnes } 9173a40ed3dSBarry Smith PetscFunctionReturn(0); 9183501a2bdSLois Curfman McInnes } 9193501a2bdSLois Curfman McInnes 920fc4dec0aSBarry Smith PetscErrorCode MatTranspose_MPIDense(Mat A,MatReuse reuse,Mat *matout) 9213501a2bdSLois Curfman McInnes { 9223501a2bdSLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense*)A->data; 9233501a2bdSLois Curfman McInnes Mat_SeqDense *Aloc = (Mat_SeqDense*)a->A->data; 9243501a2bdSLois Curfman McInnes Mat B; 925d0f46423SBarry Smith PetscInt M = A->rmap->N,N = A->cmap->N,m,n,*rwork,rstart = A->rmap->rstart; 9266849ba73SBarry Smith PetscErrorCode ierr; 92713f74950SBarry Smith PetscInt j,i; 92887828ca2SBarry Smith PetscScalar *v; 9293501a2bdSLois Curfman McInnes 9303a40ed3dSBarry Smith PetscFunctionBegin; 931cf37664fSBarry Smith if (reuse == MAT_INPLACE_MATRIX && M != N) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Supports square matrix only in-place"); 932cf37664fSBarry Smith if (reuse == MAT_INITIAL_MATRIX || reuse == MAT_INPLACE_MATRIX) { 933ce94432eSBarry Smith ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr); 934d0f46423SBarry Smith ierr = MatSetSizes(B,A->cmap->n,A->rmap->n,N,M);CHKERRQ(ierr); 9357adad957SLisandro Dalcin ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr); 9360298fd71SBarry Smith ierr = MatMPIDenseSetPreallocation(B,NULL);CHKERRQ(ierr); 937fc4dec0aSBarry Smith } else { 938fc4dec0aSBarry Smith B = *matout; 939fc4dec0aSBarry Smith } 9403501a2bdSLois Curfman McInnes 941d0f46423SBarry Smith m = a->A->rmap->n; n = a->A->cmap->n; v = Aloc->v; 942785e854fSJed Brown ierr = PetscMalloc1(m,&rwork);CHKERRQ(ierr); 9433501a2bdSLois Curfman McInnes for (i=0; i<m; i++) rwork[i] = rstart + i; 9441acff37aSSatish Balay for (j=0; j<n; j++) { 9453501a2bdSLois Curfman McInnes ierr = MatSetValues(B,1,&j,m,rwork,v,INSERT_VALUES);CHKERRQ(ierr); 9463501a2bdSLois Curfman McInnes v += m; 9473501a2bdSLois Curfman McInnes } 948606d414cSSatish Balay ierr = PetscFree(rwork);CHKERRQ(ierr); 9496d4a8577SBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 9506d4a8577SBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 951cf37664fSBarry Smith if (reuse == MAT_INITIAL_MATRIX || reuse == MAT_REUSE_MATRIX) { 9523501a2bdSLois Curfman McInnes *matout = B; 9533501a2bdSLois Curfman McInnes } else { 95428be2f97SBarry Smith ierr = MatHeaderMerge(A,&B);CHKERRQ(ierr); 9553501a2bdSLois Curfman McInnes } 9563a40ed3dSBarry Smith PetscFunctionReturn(0); 957096963f5SLois Curfman McInnes } 958096963f5SLois Curfman McInnes 95944cd7ae7SLois Curfman McInnes 9606849ba73SBarry Smith static PetscErrorCode MatDuplicate_MPIDense(Mat,MatDuplicateOption,Mat*); 961d84338a6SBarry Smith extern PetscErrorCode MatScale_MPIDense(Mat,PetscScalar); 9628965ea79SLois Curfman McInnes 9634994cf47SJed Brown PetscErrorCode MatSetUp_MPIDense(Mat A) 964273d9f13SBarry Smith { 965dfbe8321SBarry Smith PetscErrorCode ierr; 966273d9f13SBarry Smith 967273d9f13SBarry Smith PetscFunctionBegin; 968273d9f13SBarry Smith ierr = MatMPIDenseSetPreallocation(A,0);CHKERRQ(ierr); 969273d9f13SBarry Smith PetscFunctionReturn(0); 970273d9f13SBarry Smith } 971273d9f13SBarry Smith 972488007eeSBarry Smith PetscErrorCode MatAXPY_MPIDense(Mat Y,PetscScalar alpha,Mat X,MatStructure str) 973488007eeSBarry Smith { 974488007eeSBarry Smith PetscErrorCode ierr; 975488007eeSBarry Smith Mat_MPIDense *A = (Mat_MPIDense*)Y->data, *B = (Mat_MPIDense*)X->data; 976488007eeSBarry Smith 977488007eeSBarry Smith PetscFunctionBegin; 978488007eeSBarry Smith ierr = MatAXPY(A->A,alpha,B->A,str);CHKERRQ(ierr); 979a3fa217bSJose E. Roman ierr = PetscObjectStateIncrease((PetscObject)Y);CHKERRQ(ierr); 980488007eeSBarry Smith PetscFunctionReturn(0); 981488007eeSBarry Smith } 982488007eeSBarry Smith 9837087cfbeSBarry Smith PetscErrorCode MatConjugate_MPIDense(Mat mat) 984ba337c44SJed Brown { 985ba337c44SJed Brown Mat_MPIDense *a = (Mat_MPIDense*)mat->data; 986ba337c44SJed Brown PetscErrorCode ierr; 987ba337c44SJed Brown 988ba337c44SJed Brown PetscFunctionBegin; 989ba337c44SJed Brown ierr = MatConjugate(a->A);CHKERRQ(ierr); 990ba337c44SJed Brown PetscFunctionReturn(0); 991ba337c44SJed Brown } 992ba337c44SJed Brown 993ba337c44SJed Brown PetscErrorCode MatRealPart_MPIDense(Mat A) 994ba337c44SJed Brown { 995ba337c44SJed Brown Mat_MPIDense *a = (Mat_MPIDense*)A->data; 996ba337c44SJed Brown PetscErrorCode ierr; 997ba337c44SJed Brown 998ba337c44SJed Brown PetscFunctionBegin; 999ba337c44SJed Brown ierr = MatRealPart(a->A);CHKERRQ(ierr); 1000ba337c44SJed Brown PetscFunctionReturn(0); 1001ba337c44SJed Brown } 1002ba337c44SJed Brown 1003ba337c44SJed Brown PetscErrorCode MatImaginaryPart_MPIDense(Mat A) 1004ba337c44SJed Brown { 1005ba337c44SJed Brown Mat_MPIDense *a = (Mat_MPIDense*)A->data; 1006ba337c44SJed Brown PetscErrorCode ierr; 1007ba337c44SJed Brown 1008ba337c44SJed Brown PetscFunctionBegin; 1009ba337c44SJed Brown ierr = MatImaginaryPart(a->A);CHKERRQ(ierr); 1010ba337c44SJed Brown PetscFunctionReturn(0); 1011ba337c44SJed Brown } 1012ba337c44SJed Brown 10130716a85fSBarry Smith extern PetscErrorCode MatGetColumnNorms_SeqDense(Mat,NormType,PetscReal*); 10140716a85fSBarry Smith PetscErrorCode MatGetColumnNorms_MPIDense(Mat A,NormType type,PetscReal *norms) 10150716a85fSBarry Smith { 10160716a85fSBarry Smith PetscErrorCode ierr; 10170716a85fSBarry Smith PetscInt i,n; 10180716a85fSBarry Smith Mat_MPIDense *a = (Mat_MPIDense*) A->data; 10190716a85fSBarry Smith PetscReal *work; 10200716a85fSBarry Smith 10210716a85fSBarry Smith PetscFunctionBegin; 10220298fd71SBarry Smith ierr = MatGetSize(A,NULL,&n);CHKERRQ(ierr); 1023785e854fSJed Brown ierr = PetscMalloc1(n,&work);CHKERRQ(ierr); 10240716a85fSBarry Smith ierr = MatGetColumnNorms_SeqDense(a->A,type,work);CHKERRQ(ierr); 10250716a85fSBarry Smith if (type == NORM_2) { 10260716a85fSBarry Smith for (i=0; i<n; i++) work[i] *= work[i]; 10270716a85fSBarry Smith } 10280716a85fSBarry Smith if (type == NORM_INFINITY) { 1029b2566f29SBarry Smith ierr = MPIU_Allreduce(work,norms,n,MPIU_REAL,MPIU_MAX,A->hdr.comm);CHKERRQ(ierr); 10300716a85fSBarry Smith } else { 1031b2566f29SBarry Smith ierr = MPIU_Allreduce(work,norms,n,MPIU_REAL,MPIU_SUM,A->hdr.comm);CHKERRQ(ierr); 10320716a85fSBarry Smith } 10330716a85fSBarry Smith ierr = PetscFree(work);CHKERRQ(ierr); 10340716a85fSBarry Smith if (type == NORM_2) { 10358f1a2a5eSBarry Smith for (i=0; i<n; i++) norms[i] = PetscSqrtReal(norms[i]); 10360716a85fSBarry Smith } 10370716a85fSBarry Smith PetscFunctionReturn(0); 10380716a85fSBarry Smith } 10390716a85fSBarry Smith 104073a71a0fSBarry Smith static PetscErrorCode MatSetRandom_MPIDense(Mat x,PetscRandom rctx) 104173a71a0fSBarry Smith { 104273a71a0fSBarry Smith Mat_MPIDense *d = (Mat_MPIDense*)x->data; 104373a71a0fSBarry Smith PetscErrorCode ierr; 104473a71a0fSBarry Smith PetscScalar *a; 104573a71a0fSBarry Smith PetscInt m,n,i; 104673a71a0fSBarry Smith 104773a71a0fSBarry Smith PetscFunctionBegin; 104873a71a0fSBarry Smith ierr = MatGetSize(d->A,&m,&n);CHKERRQ(ierr); 10498c778c55SBarry Smith ierr = MatDenseGetArray(d->A,&a);CHKERRQ(ierr); 105073a71a0fSBarry Smith for (i=0; i<m*n; i++) { 105173a71a0fSBarry Smith ierr = PetscRandomGetValue(rctx,a+i);CHKERRQ(ierr); 105273a71a0fSBarry Smith } 10538c778c55SBarry Smith ierr = MatDenseRestoreArray(d->A,&a);CHKERRQ(ierr); 105473a71a0fSBarry Smith PetscFunctionReturn(0); 105573a71a0fSBarry Smith } 105673a71a0fSBarry Smith 1057fd4e9aacSBarry Smith extern PetscErrorCode MatMatMultNumeric_MPIDense(Mat A,Mat,Mat); 1058fd4e9aacSBarry Smith 10593b49f96aSBarry Smith static PetscErrorCode MatMissingDiagonal_MPIDense(Mat A,PetscBool *missing,PetscInt *d) 10603b49f96aSBarry Smith { 10613b49f96aSBarry Smith PetscFunctionBegin; 10623b49f96aSBarry Smith *missing = PETSC_FALSE; 10633b49f96aSBarry Smith PetscFunctionReturn(0); 10643b49f96aSBarry Smith } 10653b49f96aSBarry Smith 10668965ea79SLois Curfman McInnes /* -------------------------------------------------------------------*/ 106709dc0095SBarry Smith static struct _MatOps MatOps_Values = { MatSetValues_MPIDense, 106809dc0095SBarry Smith MatGetRow_MPIDense, 106909dc0095SBarry Smith MatRestoreRow_MPIDense, 107009dc0095SBarry Smith MatMult_MPIDense, 107197304618SKris Buschelman /* 4*/ MatMultAdd_MPIDense, 10727c922b88SBarry Smith MatMultTranspose_MPIDense, 10737c922b88SBarry Smith MatMultTransposeAdd_MPIDense, 10748965ea79SLois Curfman McInnes 0, 107509dc0095SBarry Smith 0, 107609dc0095SBarry Smith 0, 107797304618SKris Buschelman /* 10*/ 0, 107809dc0095SBarry Smith 0, 107909dc0095SBarry Smith 0, 108009dc0095SBarry Smith 0, 108109dc0095SBarry Smith MatTranspose_MPIDense, 108297304618SKris Buschelman /* 15*/ MatGetInfo_MPIDense, 10836e4ee0c6SHong Zhang MatEqual_MPIDense, 108409dc0095SBarry Smith MatGetDiagonal_MPIDense, 10855b2fa520SLois Curfman McInnes MatDiagonalScale_MPIDense, 108609dc0095SBarry Smith MatNorm_MPIDense, 108797304618SKris Buschelman /* 20*/ MatAssemblyBegin_MPIDense, 108809dc0095SBarry Smith MatAssemblyEnd_MPIDense, 108909dc0095SBarry Smith MatSetOption_MPIDense, 109009dc0095SBarry Smith MatZeroEntries_MPIDense, 1091d519adbfSMatthew Knepley /* 24*/ MatZeroRows_MPIDense, 1092919b68f7SBarry Smith 0, 109301b82886SBarry Smith 0, 109401b82886SBarry Smith 0, 109501b82886SBarry Smith 0, 10964994cf47SJed Brown /* 29*/ MatSetUp_MPIDense, 1097273d9f13SBarry Smith 0, 109809dc0095SBarry Smith 0, 1099c56a70eeSBarry Smith MatGetDiagonalBlock_MPIDense, 11008c778c55SBarry Smith 0, 1101d519adbfSMatthew Knepley /* 34*/ MatDuplicate_MPIDense, 110209dc0095SBarry Smith 0, 110309dc0095SBarry Smith 0, 110409dc0095SBarry Smith 0, 110509dc0095SBarry Smith 0, 1106d519adbfSMatthew Knepley /* 39*/ MatAXPY_MPIDense, 11077dae84e0SHong Zhang MatCreateSubMatrices_MPIDense, 110809dc0095SBarry Smith 0, 110909dc0095SBarry Smith MatGetValues_MPIDense, 111009dc0095SBarry Smith 0, 1111d519adbfSMatthew Knepley /* 44*/ 0, 111209dc0095SBarry Smith MatScale_MPIDense, 11137d68702bSBarry Smith MatShift_Basic, 111409dc0095SBarry Smith 0, 111509dc0095SBarry Smith 0, 111673a71a0fSBarry Smith /* 49*/ MatSetRandom_MPIDense, 111709dc0095SBarry Smith 0, 111809dc0095SBarry Smith 0, 111909dc0095SBarry Smith 0, 112009dc0095SBarry Smith 0, 1121d519adbfSMatthew Knepley /* 54*/ 0, 112209dc0095SBarry Smith 0, 112309dc0095SBarry Smith 0, 112409dc0095SBarry Smith 0, 112509dc0095SBarry Smith 0, 11267dae84e0SHong Zhang /* 59*/ MatCreateSubMatrix_MPIDense, 1127b9b97703SBarry Smith MatDestroy_MPIDense, 1128b9b97703SBarry Smith MatView_MPIDense, 1129357abbc8SBarry Smith 0, 113097304618SKris Buschelman 0, 1131d519adbfSMatthew Knepley /* 64*/ 0, 113297304618SKris Buschelman 0, 113397304618SKris Buschelman 0, 113497304618SKris Buschelman 0, 113597304618SKris Buschelman 0, 1136d519adbfSMatthew Knepley /* 69*/ 0, 113797304618SKris Buschelman 0, 113897304618SKris Buschelman 0, 113997304618SKris Buschelman 0, 114097304618SKris Buschelman 0, 1141d519adbfSMatthew Knepley /* 74*/ 0, 114297304618SKris Buschelman 0, 114397304618SKris Buschelman 0, 114497304618SKris Buschelman 0, 114597304618SKris Buschelman 0, 1146d519adbfSMatthew Knepley /* 79*/ 0, 114797304618SKris Buschelman 0, 114897304618SKris Buschelman 0, 114997304618SKris Buschelman 0, 11505bba2384SShri Abhyankar /* 83*/ MatLoad_MPIDense, 1151865e5f61SKris Buschelman 0, 1152865e5f61SKris Buschelman 0, 1153865e5f61SKris Buschelman 0, 1154865e5f61SKris Buschelman 0, 1155865e5f61SKris Buschelman 0, 11569775878aSHong Zhang #if defined(PETSC_HAVE_ELEMENTAL) 1157320f2790SHong Zhang /* 89*/ MatMatMult_MPIDense_MPIDense, 1158320f2790SHong Zhang MatMatMultSymbolic_MPIDense_MPIDense, 11599775878aSHong Zhang #else 11609775878aSHong Zhang /* 89*/ 0, 1161865e5f61SKris Buschelman 0, 11629775878aSHong Zhang #endif 1163fd4e9aacSBarry Smith MatMatMultNumeric_MPIDense, 11642fbe02b9SBarry Smith 0, 1165ba337c44SJed Brown 0, 1166d519adbfSMatthew Knepley /* 94*/ 0, 1167865e5f61SKris Buschelman 0, 1168865e5f61SKris Buschelman 0, 1169ba337c44SJed Brown 0, 1170ba337c44SJed Brown 0, 1171ba337c44SJed Brown /* 99*/ 0, 1172ba337c44SJed Brown 0, 1173ba337c44SJed Brown 0, 1174ba337c44SJed Brown MatConjugate_MPIDense, 1175ba337c44SJed Brown 0, 1176ba337c44SJed Brown /*104*/ 0, 1177ba337c44SJed Brown MatRealPart_MPIDense, 1178ba337c44SJed Brown MatImaginaryPart_MPIDense, 117986d161a7SShri Abhyankar 0, 118086d161a7SShri Abhyankar 0, 118186d161a7SShri Abhyankar /*109*/ 0, 118286d161a7SShri Abhyankar 0, 118386d161a7SShri Abhyankar 0, 118486d161a7SShri Abhyankar 0, 11853b49f96aSBarry Smith MatMissingDiagonal_MPIDense, 118686d161a7SShri Abhyankar /*114*/ 0, 118786d161a7SShri Abhyankar 0, 118886d161a7SShri Abhyankar 0, 118986d161a7SShri Abhyankar 0, 119086d161a7SShri Abhyankar 0, 119186d161a7SShri Abhyankar /*119*/ 0, 119286d161a7SShri Abhyankar 0, 119386d161a7SShri Abhyankar 0, 11940716a85fSBarry Smith 0, 11950716a85fSBarry Smith 0, 11960716a85fSBarry Smith /*124*/ 0, 11973964eb88SJed Brown MatGetColumnNorms_MPIDense, 11983964eb88SJed Brown 0, 11993964eb88SJed Brown 0, 12003964eb88SJed Brown 0, 12013964eb88SJed Brown /*129*/ 0, 1202cb20be35SHong Zhang MatTransposeMatMult_MPIDense_MPIDense, 1203cb20be35SHong Zhang MatTransposeMatMultSymbolic_MPIDense_MPIDense, 1204cb20be35SHong Zhang MatTransposeMatMultNumeric_MPIDense_MPIDense, 12053964eb88SJed Brown 0, 12063964eb88SJed Brown /*134*/ 0, 12073964eb88SJed Brown 0, 12083964eb88SJed Brown 0, 12093964eb88SJed Brown 0, 12103964eb88SJed Brown 0, 12113964eb88SJed Brown /*139*/ 0, 1212f9426fe0SMark Adams 0, 12133964eb88SJed Brown 0 1214ba337c44SJed Brown }; 12158965ea79SLois Curfman McInnes 12167087cfbeSBarry Smith PetscErrorCode MatMPIDenseSetPreallocation_MPIDense(Mat mat,PetscScalar *data) 1217a23d5eceSKris Buschelman { 1218a23d5eceSKris Buschelman Mat_MPIDense *a; 1219dfbe8321SBarry Smith PetscErrorCode ierr; 1220a23d5eceSKris Buschelman 1221a23d5eceSKris Buschelman PetscFunctionBegin; 1222a23d5eceSKris Buschelman mat->preallocated = PETSC_TRUE; 1223a23d5eceSKris Buschelman /* Note: For now, when data is specified above, this assumes the user correctly 1224a23d5eceSKris Buschelman allocates the local dense storage space. We should add error checking. */ 1225a23d5eceSKris Buschelman 1226a23d5eceSKris Buschelman a = (Mat_MPIDense*)mat->data; 122734ef9618SShri Abhyankar ierr = PetscLayoutSetUp(mat->rmap);CHKERRQ(ierr); 122834ef9618SShri Abhyankar ierr = PetscLayoutSetUp(mat->cmap);CHKERRQ(ierr); 122934ef9618SShri Abhyankar a->nvec = mat->cmap->n; 123034ef9618SShri Abhyankar 1231f69a0ea3SMatthew Knepley ierr = MatCreate(PETSC_COMM_SELF,&a->A);CHKERRQ(ierr); 1232d0f46423SBarry Smith ierr = MatSetSizes(a->A,mat->rmap->n,mat->cmap->N,mat->rmap->n,mat->cmap->N);CHKERRQ(ierr); 12335c5985e7SKris Buschelman ierr = MatSetType(a->A,MATSEQDENSE);CHKERRQ(ierr); 12345c5985e7SKris Buschelman ierr = MatSeqDenseSetPreallocation(a->A,data);CHKERRQ(ierr); 12353bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)a->A);CHKERRQ(ierr); 1236a23d5eceSKris Buschelman PetscFunctionReturn(0); 1237a23d5eceSKris Buschelman } 1238a23d5eceSKris Buschelman 123965b80a83SHong Zhang #if defined(PETSC_HAVE_ELEMENTAL) 1240cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatConvert_MPIDense_Elemental(Mat A, MatType newtype,MatReuse reuse,Mat *newmat) 12418baccfbdSHong Zhang { 12428ea901baSHong Zhang Mat mat_elemental; 12438ea901baSHong Zhang PetscErrorCode ierr; 124432d7a744SHong Zhang PetscScalar *v; 124532d7a744SHong Zhang PetscInt m=A->rmap->n,N=A->cmap->N,rstart=A->rmap->rstart,i,*rows,*cols; 12468ea901baSHong Zhang 12478baccfbdSHong Zhang PetscFunctionBegin; 1248378336b6SHong Zhang if (reuse == MAT_REUSE_MATRIX) { 1249378336b6SHong Zhang mat_elemental = *newmat; 1250378336b6SHong Zhang ierr = MatZeroEntries(*newmat);CHKERRQ(ierr); 1251378336b6SHong Zhang } else { 1252378336b6SHong Zhang ierr = MatCreate(PetscObjectComm((PetscObject)A), &mat_elemental);CHKERRQ(ierr); 1253378336b6SHong Zhang ierr = MatSetSizes(mat_elemental,PETSC_DECIDE,PETSC_DECIDE,A->rmap->N,A->cmap->N);CHKERRQ(ierr); 1254378336b6SHong Zhang ierr = MatSetType(mat_elemental,MATELEMENTAL);CHKERRQ(ierr); 1255378336b6SHong Zhang ierr = MatSetUp(mat_elemental);CHKERRQ(ierr); 125632d7a744SHong Zhang ierr = MatSetOption(mat_elemental,MAT_ROW_ORIENTED,PETSC_FALSE);CHKERRQ(ierr); 1257378336b6SHong Zhang } 1258378336b6SHong Zhang 125932d7a744SHong Zhang ierr = PetscMalloc2(m,&rows,N,&cols);CHKERRQ(ierr); 126032d7a744SHong Zhang for (i=0; i<N; i++) cols[i] = i; 126132d7a744SHong Zhang for (i=0; i<m; i++) rows[i] = rstart + i; 12628ea901baSHong Zhang 12638ea901baSHong Zhang /* PETSc-Elemental interaface uses axpy for setting off-processor entries, only ADD_VALUES is allowed */ 126432d7a744SHong Zhang ierr = MatDenseGetArray(A,&v);CHKERRQ(ierr); 126532d7a744SHong Zhang ierr = MatSetValues(mat_elemental,m,rows,N,cols,v,ADD_VALUES);CHKERRQ(ierr); 12668ea901baSHong Zhang ierr = MatAssemblyBegin(mat_elemental, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 12678ea901baSHong Zhang ierr = MatAssemblyEnd(mat_elemental, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 126832d7a744SHong Zhang ierr = MatDenseRestoreArray(A,&v);CHKERRQ(ierr); 126932d7a744SHong Zhang ierr = PetscFree2(rows,cols);CHKERRQ(ierr); 12708ea901baSHong Zhang 1271511c6705SHong Zhang if (reuse == MAT_INPLACE_MATRIX) { 127228be2f97SBarry Smith ierr = MatHeaderReplace(A,&mat_elemental);CHKERRQ(ierr); 12738ea901baSHong Zhang } else { 12748ea901baSHong Zhang *newmat = mat_elemental; 12758ea901baSHong Zhang } 12768baccfbdSHong Zhang PetscFunctionReturn(0); 12778baccfbdSHong Zhang } 127865b80a83SHong Zhang #endif 12798baccfbdSHong Zhang 1280*86aefd0dSHong Zhang static PetscErrorCode MatGetColumn_MPIDense(Mat A,PetscInt col,PetscScalar **vals) 1281*86aefd0dSHong Zhang { 1282*86aefd0dSHong Zhang Mat_MPIDense *mat = (Mat_MPIDense*)A->data; 1283*86aefd0dSHong Zhang PetscErrorCode ierr; 1284*86aefd0dSHong Zhang 1285*86aefd0dSHong Zhang PetscFunctionBegin; 1286*86aefd0dSHong Zhang ierr = MatDenseGetColumn(mat->A,col,vals);CHKERRQ(ierr); 1287*86aefd0dSHong Zhang PetscFunctionReturn(0); 1288*86aefd0dSHong Zhang } 1289*86aefd0dSHong Zhang 1290*86aefd0dSHong Zhang static PetscErrorCode MatRestoreColumn_MPIDense(Mat A,PetscScalar **vals) 1291*86aefd0dSHong Zhang { 1292*86aefd0dSHong Zhang Mat_MPIDense *mat = (Mat_MPIDense*)A->data; 1293*86aefd0dSHong Zhang PetscErrorCode ierr; 1294*86aefd0dSHong Zhang 1295*86aefd0dSHong Zhang PetscFunctionBegin; 1296*86aefd0dSHong Zhang ierr = MatDenseRestoreColumn(mat->A,vals);CHKERRQ(ierr); 1297*86aefd0dSHong Zhang PetscFunctionReturn(0); 1298*86aefd0dSHong Zhang } 1299*86aefd0dSHong Zhang 13008cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_MPIDense(Mat mat) 1301273d9f13SBarry Smith { 1302273d9f13SBarry Smith Mat_MPIDense *a; 1303dfbe8321SBarry Smith PetscErrorCode ierr; 1304273d9f13SBarry Smith 1305273d9f13SBarry Smith PetscFunctionBegin; 1306b00a9115SJed Brown ierr = PetscNewLog(mat,&a);CHKERRQ(ierr); 1307b0a32e0cSBarry Smith mat->data = (void*)a; 1308273d9f13SBarry Smith ierr = PetscMemcpy(mat->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 1309273d9f13SBarry Smith 1310273d9f13SBarry Smith mat->insertmode = NOT_SET_VALUES; 1311ce94432eSBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)mat),&a->rank);CHKERRQ(ierr); 1312ce94432eSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)mat),&a->size);CHKERRQ(ierr); 1313273d9f13SBarry Smith 1314273d9f13SBarry Smith /* build cache for off array entries formed */ 1315273d9f13SBarry Smith a->donotstash = PETSC_FALSE; 13162205254eSKarl Rupp 1317ce94432eSBarry Smith ierr = MatStashCreate_Private(PetscObjectComm((PetscObject)mat),1,&mat->stash);CHKERRQ(ierr); 1318273d9f13SBarry Smith 1319273d9f13SBarry Smith /* stuff used for matrix vector multiply */ 1320273d9f13SBarry Smith a->lvec = 0; 1321273d9f13SBarry Smith a->Mvctx = 0; 1322273d9f13SBarry Smith a->roworiented = PETSC_TRUE; 1323273d9f13SBarry Smith 1324bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatDenseGetArray_C",MatDenseGetArray_MPIDense);CHKERRQ(ierr); 1325d3042a70SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatDensePlaceArray_C",MatDensePlaceArray_MPIDense);CHKERRQ(ierr); 1326d3042a70SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatDenseResetArray_C",MatDenseResetArray_MPIDense);CHKERRQ(ierr); 1327bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatDenseRestoreArray_C",MatDenseRestoreArray_MPIDense);CHKERRQ(ierr); 13288baccfbdSHong Zhang #if defined(PETSC_HAVE_ELEMENTAL) 13298baccfbdSHong Zhang ierr = PetscObjectComposeFunction((PetscObject)mat,"MatConvert_mpidense_elemental_C",MatConvert_MPIDense_Elemental);CHKERRQ(ierr); 13308baccfbdSHong Zhang #endif 1331bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMPIDenseSetPreallocation_C",MatMPIDenseSetPreallocation_MPIDense);CHKERRQ(ierr); 1332bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMatMult_mpiaij_mpidense_C",MatMatMult_MPIAIJ_MPIDense);CHKERRQ(ierr); 1333bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMatMultSymbolic_mpiaij_mpidense_C",MatMatMultSymbolic_MPIAIJ_MPIDense);CHKERRQ(ierr); 1334bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMatMultNumeric_mpiaij_mpidense_C",MatMatMultNumeric_MPIAIJ_MPIDense);CHKERRQ(ierr); 13358949adfdSHong Zhang 13368949adfdSHong Zhang ierr = PetscObjectComposeFunction((PetscObject)mat,"MatTransposeMatMult_mpiaij_mpidense_C",MatTransposeMatMult_MPIAIJ_MPIDense);CHKERRQ(ierr); 13378949adfdSHong Zhang ierr = PetscObjectComposeFunction((PetscObject)mat,"MatTransposeMatMultSymbolic_mpiaij_mpidense_C",MatTransposeMatMultSymbolic_MPIAIJ_MPIDense);CHKERRQ(ierr); 13388949adfdSHong Zhang ierr = PetscObjectComposeFunction((PetscObject)mat,"MatTransposeMatMultNumeric_mpiaij_mpidense_C",MatTransposeMatMultNumeric_MPIAIJ_MPIDense);CHKERRQ(ierr); 1339*86aefd0dSHong Zhang ierr = PetscObjectComposeFunction((PetscObject)mat,"MatDenseGetColumn_C",MatGetColumn_MPIDense);CHKERRQ(ierr); 1340*86aefd0dSHong Zhang ierr = PetscObjectComposeFunction((PetscObject)mat,"MatDenseRestoreColumn_C",MatRestoreColumn_MPIDense);CHKERRQ(ierr); 134138aed534SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)mat,MATMPIDENSE);CHKERRQ(ierr); 1342273d9f13SBarry Smith PetscFunctionReturn(0); 1343273d9f13SBarry Smith } 1344273d9f13SBarry Smith 1345209238afSKris Buschelman /*MC 1346002d173eSKris Buschelman MATDENSE - MATDENSE = "dense" - A matrix type to be used for dense matrices. 1347209238afSKris Buschelman 1348209238afSKris Buschelman This matrix type is identical to MATSEQDENSE when constructed with a single process communicator, 1349209238afSKris Buschelman and MATMPIDENSE otherwise. 1350209238afSKris Buschelman 1351209238afSKris Buschelman Options Database Keys: 1352209238afSKris Buschelman . -mat_type dense - sets the matrix type to "dense" during a call to MatSetFromOptions() 1353209238afSKris Buschelman 1354209238afSKris Buschelman Level: beginner 1355209238afSKris Buschelman 135601b82886SBarry Smith 1357209238afSKris Buschelman .seealso: MatCreateMPIDense,MATSEQDENSE,MATMPIDENSE 1358209238afSKris Buschelman M*/ 1359209238afSKris Buschelman 1360273d9f13SBarry Smith /*@C 1361273d9f13SBarry Smith MatMPIDenseSetPreallocation - Sets the array used to store the matrix entries 1362273d9f13SBarry Smith 1363273d9f13SBarry Smith Not collective 1364273d9f13SBarry Smith 1365273d9f13SBarry Smith Input Parameters: 13661c4f3114SJed Brown . B - the matrix 13670298fd71SBarry Smith - data - optional location of matrix data. Set data=NULL for PETSc 1368273d9f13SBarry Smith to control all matrix memory allocation. 1369273d9f13SBarry Smith 1370273d9f13SBarry Smith Notes: 1371273d9f13SBarry Smith The dense format is fully compatible with standard Fortran 77 1372273d9f13SBarry Smith storage by columns. 1373273d9f13SBarry Smith 1374273d9f13SBarry Smith The data input variable is intended primarily for Fortran programmers 1375273d9f13SBarry Smith who wish to allocate their own matrix memory space. Most users should 13760298fd71SBarry Smith set data=NULL. 1377273d9f13SBarry Smith 1378273d9f13SBarry Smith Level: intermediate 1379273d9f13SBarry Smith 1380273d9f13SBarry Smith .keywords: matrix,dense, parallel 1381273d9f13SBarry Smith 1382273d9f13SBarry Smith .seealso: MatCreate(), MatCreateSeqDense(), MatSetValues() 1383273d9f13SBarry Smith @*/ 13841c4f3114SJed Brown PetscErrorCode MatMPIDenseSetPreallocation(Mat B,PetscScalar *data) 1385273d9f13SBarry Smith { 13864ac538c5SBarry Smith PetscErrorCode ierr; 1387273d9f13SBarry Smith 1388273d9f13SBarry Smith PetscFunctionBegin; 13891c4f3114SJed Brown ierr = PetscTryMethod(B,"MatMPIDenseSetPreallocation_C",(Mat,PetscScalar*),(B,data));CHKERRQ(ierr); 1390273d9f13SBarry Smith PetscFunctionReturn(0); 1391273d9f13SBarry Smith } 1392273d9f13SBarry Smith 1393d3042a70SBarry Smith /*@ 1394d3042a70SBarry Smith MatDensePlaceArray - Allows one to replace the array in a dense array with an 1395d3042a70SBarry Smith array provided by the user. This is useful to avoid copying an array 1396d3042a70SBarry Smith into a matrix 1397d3042a70SBarry Smith 1398d3042a70SBarry Smith Not Collective 1399d3042a70SBarry Smith 1400d3042a70SBarry Smith Input Parameters: 1401d3042a70SBarry Smith + mat - the matrix 1402d3042a70SBarry Smith - array - the array in column major order 1403d3042a70SBarry Smith 1404d3042a70SBarry Smith Notes: 1405d3042a70SBarry Smith You can return to the original array with a call to MatDenseResetArray(). The user is responsible for freeing this array; it will not be 1406d3042a70SBarry Smith freed when the matrix is destroyed. 1407d3042a70SBarry Smith 1408d3042a70SBarry Smith Level: developer 1409d3042a70SBarry Smith 1410d3042a70SBarry Smith .seealso: MatDenseGetArray(), MatDenseResetArray(), VecPlaceArray(), VecGetArray(), VecRestoreArray(), VecReplaceArray(), VecResetArray() 1411d3042a70SBarry Smith 1412d3042a70SBarry Smith @*/ 1413d3042a70SBarry Smith PetscErrorCode MatDensePlaceArray(Mat mat,const PetscScalar array[]) 1414d3042a70SBarry Smith { 1415d3042a70SBarry Smith PetscErrorCode ierr; 1416d3042a70SBarry Smith PetscFunctionBegin; 1417d3042a70SBarry Smith ierr = PetscUseMethod(mat,"MatDensePlaceArray_C",(Mat,const PetscScalar*),(mat,array));CHKERRQ(ierr); 1418d3042a70SBarry Smith ierr = PetscObjectStateIncrease((PetscObject)mat);CHKERRQ(ierr); 1419d3042a70SBarry Smith PetscFunctionReturn(0); 1420d3042a70SBarry Smith } 1421d3042a70SBarry Smith 1422d3042a70SBarry Smith /*@ 1423d3042a70SBarry Smith MatDenseResetArray - Resets the matrix array to that it previously had before the call to MatDensePlaceArray() 1424d3042a70SBarry Smith 1425d3042a70SBarry Smith Not Collective 1426d3042a70SBarry Smith 1427d3042a70SBarry Smith Input Parameters: 1428d3042a70SBarry Smith . mat - the matrix 1429d3042a70SBarry Smith 1430d3042a70SBarry Smith Notes: 1431d3042a70SBarry Smith You can only call this after a call to MatDensePlaceArray() 1432d3042a70SBarry Smith 1433d3042a70SBarry Smith Level: developer 1434d3042a70SBarry Smith 1435d3042a70SBarry Smith .seealso: MatDenseGetArray(), MatDensePlaceArray(), VecPlaceArray(), VecGetArray(), VecRestoreArray(), VecReplaceArray(), VecResetArray() 1436d3042a70SBarry Smith 1437d3042a70SBarry Smith @*/ 1438d3042a70SBarry Smith PetscErrorCode MatDenseResetArray(Mat mat) 1439d3042a70SBarry Smith { 1440d3042a70SBarry Smith PetscErrorCode ierr; 1441d3042a70SBarry Smith PetscFunctionBegin; 1442d3042a70SBarry Smith ierr = PetscUseMethod(mat,"MatDenseResetArray_C",(Mat),(mat));CHKERRQ(ierr); 1443d3042a70SBarry Smith ierr = PetscObjectStateIncrease((PetscObject)mat);CHKERRQ(ierr); 1444d3042a70SBarry Smith PetscFunctionReturn(0); 1445d3042a70SBarry Smith } 1446d3042a70SBarry Smith 14478965ea79SLois Curfman McInnes /*@C 144869b1f4b7SBarry Smith MatCreateDense - Creates a parallel matrix in dense format. 14498965ea79SLois Curfman McInnes 1450db81eaa0SLois Curfman McInnes Collective on MPI_Comm 1451db81eaa0SLois Curfman McInnes 14528965ea79SLois Curfman McInnes Input Parameters: 1453db81eaa0SLois Curfman McInnes + comm - MPI communicator 14548965ea79SLois Curfman McInnes . m - number of local rows (or PETSC_DECIDE to have calculated if M is given) 1455db81eaa0SLois Curfman McInnes . n - number of local columns (or PETSC_DECIDE to have calculated if N is given) 14568965ea79SLois Curfman McInnes . M - number of global rows (or PETSC_DECIDE to have calculated if m is given) 1457db81eaa0SLois Curfman McInnes . N - number of global columns (or PETSC_DECIDE to have calculated if n is given) 14586cfe35ddSJose E. Roman - data - optional location of matrix data. Set data=NULL (PETSC_NULL_SCALAR for Fortran users) for PETSc 1459dfc5480cSLois Curfman McInnes to control all matrix memory allocation. 14608965ea79SLois Curfman McInnes 14618965ea79SLois Curfman McInnes Output Parameter: 1462477f1c0bSLois Curfman McInnes . A - the matrix 14638965ea79SLois Curfman McInnes 1464b259b22eSLois Curfman McInnes Notes: 146539ddd567SLois Curfman McInnes The dense format is fully compatible with standard Fortran 77 146639ddd567SLois Curfman McInnes storage by columns. 14678965ea79SLois Curfman McInnes 146818f449edSLois Curfman McInnes The data input variable is intended primarily for Fortran programmers 146918f449edSLois Curfman McInnes who wish to allocate their own matrix memory space. Most users should 14706cfe35ddSJose E. Roman set data=NULL (PETSC_NULL_SCALAR for Fortran users). 147118f449edSLois Curfman McInnes 14728965ea79SLois Curfman McInnes The user MUST specify either the local or global matrix dimensions 14738965ea79SLois Curfman McInnes (possibly both). 14748965ea79SLois Curfman McInnes 1475027ccd11SLois Curfman McInnes Level: intermediate 1476027ccd11SLois Curfman McInnes 147739ddd567SLois Curfman McInnes .keywords: matrix,dense, parallel 14788965ea79SLois Curfman McInnes 147939ddd567SLois Curfman McInnes .seealso: MatCreate(), MatCreateSeqDense(), MatSetValues() 14808965ea79SLois Curfman McInnes @*/ 148169b1f4b7SBarry Smith PetscErrorCode MatCreateDense(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscScalar *data,Mat *A) 14828965ea79SLois Curfman McInnes { 14836849ba73SBarry Smith PetscErrorCode ierr; 148413f74950SBarry Smith PetscMPIInt size; 14858965ea79SLois Curfman McInnes 14863a40ed3dSBarry Smith PetscFunctionBegin; 1487f69a0ea3SMatthew Knepley ierr = MatCreate(comm,A);CHKERRQ(ierr); 1488f69a0ea3SMatthew Knepley ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr); 1489273d9f13SBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1490273d9f13SBarry Smith if (size > 1) { 1491273d9f13SBarry Smith ierr = MatSetType(*A,MATMPIDENSE);CHKERRQ(ierr); 1492273d9f13SBarry Smith ierr = MatMPIDenseSetPreallocation(*A,data);CHKERRQ(ierr); 14936cfe35ddSJose E. Roman if (data) { /* user provided data array, so no need to assemble */ 14946cfe35ddSJose E. Roman ierr = MatSetUpMultiply_MPIDense(*A);CHKERRQ(ierr); 14956cfe35ddSJose E. Roman (*A)->assembled = PETSC_TRUE; 14966cfe35ddSJose E. Roman } 1497273d9f13SBarry Smith } else { 1498273d9f13SBarry Smith ierr = MatSetType(*A,MATSEQDENSE);CHKERRQ(ierr); 1499273d9f13SBarry Smith ierr = MatSeqDenseSetPreallocation(*A,data);CHKERRQ(ierr); 15008c469469SLois Curfman McInnes } 15013a40ed3dSBarry Smith PetscFunctionReturn(0); 15028965ea79SLois Curfman McInnes } 15038965ea79SLois Curfman McInnes 15046849ba73SBarry Smith static PetscErrorCode MatDuplicate_MPIDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat) 15058965ea79SLois Curfman McInnes { 15068965ea79SLois Curfman McInnes Mat mat; 15073501a2bdSLois Curfman McInnes Mat_MPIDense *a,*oldmat = (Mat_MPIDense*)A->data; 1508dfbe8321SBarry Smith PetscErrorCode ierr; 15098965ea79SLois Curfman McInnes 15103a40ed3dSBarry Smith PetscFunctionBegin; 15118965ea79SLois Curfman McInnes *newmat = 0; 1512ce94432eSBarry Smith ierr = MatCreate(PetscObjectComm((PetscObject)A),&mat);CHKERRQ(ierr); 1513d0f46423SBarry Smith ierr = MatSetSizes(mat,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr); 15147adad957SLisandro Dalcin ierr = MatSetType(mat,((PetscObject)A)->type_name);CHKERRQ(ierr); 1515834f8fabSBarry Smith a = (Mat_MPIDense*)mat->data; 1516e04c1aa4SHong Zhang ierr = PetscMemcpy(mat->ops,A->ops,sizeof(struct _MatOps));CHKERRQ(ierr); 15175aa7edbeSHong Zhang 1518d5f3da31SBarry Smith mat->factortype = A->factortype; 1519c456f294SBarry Smith mat->assembled = PETSC_TRUE; 1520273d9f13SBarry Smith mat->preallocated = PETSC_TRUE; 15218965ea79SLois Curfman McInnes 15228965ea79SLois Curfman McInnes a->size = oldmat->size; 15238965ea79SLois Curfman McInnes a->rank = oldmat->rank; 1524e0fa3b82SLois Curfman McInnes mat->insertmode = NOT_SET_VALUES; 1525b9b97703SBarry Smith a->nvec = oldmat->nvec; 15263782ba37SSatish Balay a->donotstash = oldmat->donotstash; 1527e04c1aa4SHong Zhang 15281e1e43feSBarry Smith ierr = PetscLayoutReference(A->rmap,&mat->rmap);CHKERRQ(ierr); 15291e1e43feSBarry Smith ierr = PetscLayoutReference(A->cmap,&mat->cmap);CHKERRQ(ierr); 15308965ea79SLois Curfman McInnes 1531329f5518SBarry Smith ierr = MatSetUpMultiply_MPIDense(mat);CHKERRQ(ierr); 15325609ef8eSBarry Smith ierr = MatDuplicate(oldmat->A,cpvalues,&a->A);CHKERRQ(ierr); 15333bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)a->A);CHKERRQ(ierr); 153401b82886SBarry Smith 15358965ea79SLois Curfman McInnes *newmat = mat; 15363a40ed3dSBarry Smith PetscFunctionReturn(0); 15378965ea79SLois Curfman McInnes } 15388965ea79SLois Curfman McInnes 15399d36ed5fSBarry Smith PetscErrorCode MatLoad_MPIDense_DenseInFile(MPI_Comm comm,PetscInt fd,PetscInt M,PetscInt N,Mat newmat) 154086d161a7SShri Abhyankar { 154186d161a7SShri Abhyankar PetscErrorCode ierr; 154286d161a7SShri Abhyankar PetscMPIInt rank,size; 154374c13d95SBarry Smith const PetscInt *rowners; 154474c13d95SBarry Smith PetscInt i,m,n,nz,j,mMax; 154586d161a7SShri Abhyankar PetscScalar *array,*vals,*vals_ptr; 15469d36ed5fSBarry Smith Mat_MPIDense *a = (Mat_MPIDense*)newmat->data; 154786d161a7SShri Abhyankar 154886d161a7SShri Abhyankar PetscFunctionBegin; 154986d161a7SShri Abhyankar ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 155086d161a7SShri Abhyankar ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 155186d161a7SShri Abhyankar 155274c13d95SBarry Smith /* determine ownership of rows and columns */ 155374c13d95SBarry Smith m = (newmat->rmap->n < 0) ? PETSC_DECIDE : newmat->rmap->n; 155474c13d95SBarry Smith n = (newmat->cmap->n < 0) ? PETSC_DECIDE : newmat->cmap->n; 155586d161a7SShri Abhyankar 155674c13d95SBarry Smith ierr = MatSetSizes(newmat,m,n,M,N);CHKERRQ(ierr); 15579d36ed5fSBarry Smith if (!a->A || !((Mat_SeqDense*)(a->A->data))->user_alloc) { 15580298fd71SBarry Smith ierr = MatMPIDenseSetPreallocation(newmat,NULL);CHKERRQ(ierr); 15599d36ed5fSBarry Smith } 15608c778c55SBarry Smith ierr = MatDenseGetArray(newmat,&array);CHKERRQ(ierr); 156149467e7dSBarry Smith ierr = MatGetLocalSize(newmat,&m,NULL);CHKERRQ(ierr); 156274c13d95SBarry Smith ierr = MatGetOwnershipRanges(newmat,&rowners);CHKERRQ(ierr); 1563396e826eSBarry Smith ierr = MPI_Reduce(&m,&mMax,1,MPIU_INT,MPI_MAX,0,comm);CHKERRQ(ierr); 156486d161a7SShri Abhyankar if (!rank) { 15659d36ed5fSBarry Smith ierr = PetscMalloc1(mMax*N,&vals);CHKERRQ(ierr); 156686d161a7SShri Abhyankar 156786d161a7SShri Abhyankar /* read in my part of the matrix numerical values */ 156886d161a7SShri Abhyankar ierr = PetscBinaryRead(fd,vals,m*N,PETSC_SCALAR);CHKERRQ(ierr); 156986d161a7SShri Abhyankar 157086d161a7SShri Abhyankar /* insert into matrix-by row (this is why cannot directly read into array */ 157186d161a7SShri Abhyankar vals_ptr = vals; 157286d161a7SShri Abhyankar for (i=0; i<m; i++) { 157386d161a7SShri Abhyankar for (j=0; j<N; j++) { 157486d161a7SShri Abhyankar array[i + j*m] = *vals_ptr++; 157586d161a7SShri Abhyankar } 157686d161a7SShri Abhyankar } 157786d161a7SShri Abhyankar 157886d161a7SShri Abhyankar /* read in other processors and ship out */ 157986d161a7SShri Abhyankar for (i=1; i<size; i++) { 158086d161a7SShri Abhyankar nz = (rowners[i+1] - rowners[i])*N; 158186d161a7SShri Abhyankar ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr); 1582a25532f0SBarry Smith ierr = MPIULong_Send(vals,nz,MPIU_SCALAR,i,((PetscObject)(newmat))->tag,comm);CHKERRQ(ierr); 158386d161a7SShri Abhyankar } 158486d161a7SShri Abhyankar } else { 158586d161a7SShri Abhyankar /* receive numeric values */ 1586785e854fSJed Brown ierr = PetscMalloc1(m*N,&vals);CHKERRQ(ierr); 158786d161a7SShri Abhyankar 158886d161a7SShri Abhyankar /* receive message of values*/ 1589a25532f0SBarry Smith ierr = MPIULong_Recv(vals,m*N,MPIU_SCALAR,0,((PetscObject)(newmat))->tag,comm);CHKERRQ(ierr); 159086d161a7SShri Abhyankar 159186d161a7SShri Abhyankar /* insert into matrix-by row (this is why cannot directly read into array */ 159286d161a7SShri Abhyankar vals_ptr = vals; 159386d161a7SShri Abhyankar for (i=0; i<m; i++) { 159486d161a7SShri Abhyankar for (j=0; j<N; j++) { 159586d161a7SShri Abhyankar array[i + j*m] = *vals_ptr++; 159686d161a7SShri Abhyankar } 159786d161a7SShri Abhyankar } 159886d161a7SShri Abhyankar } 15998c778c55SBarry Smith ierr = MatDenseRestoreArray(newmat,&array);CHKERRQ(ierr); 160086d161a7SShri Abhyankar ierr = PetscFree(vals);CHKERRQ(ierr); 160186d161a7SShri Abhyankar ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 160286d161a7SShri Abhyankar ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 160386d161a7SShri Abhyankar PetscFunctionReturn(0); 160486d161a7SShri Abhyankar } 160586d161a7SShri Abhyankar 1606112444f4SShri Abhyankar PetscErrorCode MatLoad_MPIDense(Mat newmat,PetscViewer viewer) 160786d161a7SShri Abhyankar { 1608dfdea288SBarry Smith Mat_MPIDense *a; 160986d161a7SShri Abhyankar PetscScalar *vals,*svals; 1610ce94432eSBarry Smith MPI_Comm comm; 161186d161a7SShri Abhyankar MPI_Status status; 161249467e7dSBarry Smith PetscMPIInt rank,size,tag = ((PetscObject)viewer)->tag,*rowners,*sndcounts,m,n,maxnz; 161386d161a7SShri Abhyankar PetscInt header[4],*rowlengths = 0,M,N,*cols; 161449467e7dSBarry Smith PetscInt *ourlens,*procsnz = 0,jj,*mycols,*smycols; 16159d36ed5fSBarry Smith PetscInt i,nz,j,rstart,rend; 161686d161a7SShri Abhyankar int fd; 161786d161a7SShri Abhyankar PetscErrorCode ierr; 161886d161a7SShri Abhyankar 161986d161a7SShri Abhyankar PetscFunctionBegin; 1620c98fd787SBarry Smith /* force binary viewer to load .info file if it has not yet done so */ 1621c98fd787SBarry Smith ierr = PetscViewerSetUp(viewer);CHKERRQ(ierr); 1622ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)viewer,&comm);CHKERRQ(ierr); 162386d161a7SShri Abhyankar ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 162486d161a7SShri Abhyankar ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 162586d161a7SShri Abhyankar ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 16265872f025SBarry Smith if (!rank) { 162786d161a7SShri Abhyankar ierr = PetscBinaryRead(fd,(char*)header,4,PETSC_INT);CHKERRQ(ierr); 162886d161a7SShri Abhyankar if (header[0] != MAT_FILE_CLASSID) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"not matrix object"); 162986d161a7SShri Abhyankar } 163086d161a7SShri Abhyankar ierr = MPI_Bcast(header+1,3,MPIU_INT,0,comm);CHKERRQ(ierr); 163186d161a7SShri Abhyankar M = header[1]; N = header[2]; nz = header[3]; 163286d161a7SShri Abhyankar 163386d161a7SShri Abhyankar /* If global rows/cols are set to PETSC_DECIDE, set it to the sizes given in the file */ 16349d36ed5fSBarry Smith if (newmat->rmap->N < 0) newmat->rmap->N = M; 16359d36ed5fSBarry Smith if (newmat->cmap->N < 0) newmat->cmap->N = N; 163686d161a7SShri Abhyankar 16379d36ed5fSBarry Smith if (newmat->rmap->N != M) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED, "Inconsistent # of rows:Matrix in file has (%D) and input matrix has (%D)",M,newmat->rmap->N); 16389d36ed5fSBarry Smith if (newmat->cmap->N != N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED, "Inconsistent # of cols:Matrix in file has (%D) and input matrix has (%D)",N,newmat->cmap->N); 163986d161a7SShri Abhyankar 164086d161a7SShri Abhyankar /* 164186d161a7SShri Abhyankar Handle case where matrix is stored on disk as a dense matrix 164286d161a7SShri Abhyankar */ 164386d161a7SShri Abhyankar if (nz == MATRIX_BINARY_FORMAT_DENSE) { 16449d36ed5fSBarry Smith ierr = MatLoad_MPIDense_DenseInFile(comm,fd,M,N,newmat);CHKERRQ(ierr); 164586d161a7SShri Abhyankar PetscFunctionReturn(0); 164686d161a7SShri Abhyankar } 164786d161a7SShri Abhyankar 164886d161a7SShri Abhyankar /* determine ownership of all rows */ 16492205254eSKarl Rupp if (newmat->rmap->n < 0) { 16502205254eSKarl Rupp ierr = PetscMPIIntCast(M/size + ((M % size) > rank),&m);CHKERRQ(ierr); 16512205254eSKarl Rupp } else { 16522205254eSKarl Rupp ierr = PetscMPIIntCast(newmat->rmap->n,&m);CHKERRQ(ierr); 16532205254eSKarl Rupp } 165449467e7dSBarry Smith if (newmat->cmap->n < 0) { 165549467e7dSBarry Smith n = PETSC_DECIDE; 165649467e7dSBarry Smith } else { 165749467e7dSBarry Smith ierr = PetscMPIIntCast(newmat->cmap->n,&n);CHKERRQ(ierr); 165849467e7dSBarry Smith } 165949467e7dSBarry Smith 1660854ce69bSBarry Smith ierr = PetscMalloc1(size+2,&rowners);CHKERRQ(ierr); 166186d161a7SShri Abhyankar ierr = MPI_Allgather(&m,1,MPI_INT,rowners+1,1,MPI_INT,comm);CHKERRQ(ierr); 166286d161a7SShri Abhyankar rowners[0] = 0; 166386d161a7SShri Abhyankar for (i=2; i<=size; i++) { 166486d161a7SShri Abhyankar rowners[i] += rowners[i-1]; 166586d161a7SShri Abhyankar } 166686d161a7SShri Abhyankar rstart = rowners[rank]; 166786d161a7SShri Abhyankar rend = rowners[rank+1]; 166886d161a7SShri Abhyankar 166986d161a7SShri Abhyankar /* distribute row lengths to all processors */ 167049467e7dSBarry Smith ierr = PetscMalloc1(rend-rstart,&ourlens);CHKERRQ(ierr); 167186d161a7SShri Abhyankar if (!rank) { 1672785e854fSJed Brown ierr = PetscMalloc1(M,&rowlengths);CHKERRQ(ierr); 167386d161a7SShri Abhyankar ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr); 1674785e854fSJed Brown ierr = PetscMalloc1(size,&sndcounts);CHKERRQ(ierr); 167586d161a7SShri Abhyankar for (i=0; i<size; i++) sndcounts[i] = rowners[i+1] - rowners[i]; 167686d161a7SShri Abhyankar ierr = MPI_Scatterv(rowlengths,sndcounts,rowners,MPIU_INT,ourlens,rend-rstart,MPIU_INT,0,comm);CHKERRQ(ierr); 167786d161a7SShri Abhyankar ierr = PetscFree(sndcounts);CHKERRQ(ierr); 167886d161a7SShri Abhyankar } else { 167986d161a7SShri Abhyankar ierr = MPI_Scatterv(0,0,0,MPIU_INT,ourlens,rend-rstart,MPIU_INT,0,comm);CHKERRQ(ierr); 168086d161a7SShri Abhyankar } 168186d161a7SShri Abhyankar 168286d161a7SShri Abhyankar if (!rank) { 168386d161a7SShri Abhyankar /* calculate the number of nonzeros on each processor */ 1684785e854fSJed Brown ierr = PetscMalloc1(size,&procsnz);CHKERRQ(ierr); 168586d161a7SShri Abhyankar ierr = PetscMemzero(procsnz,size*sizeof(PetscInt));CHKERRQ(ierr); 168686d161a7SShri Abhyankar for (i=0; i<size; i++) { 168786d161a7SShri Abhyankar for (j=rowners[i]; j< rowners[i+1]; j++) { 168886d161a7SShri Abhyankar procsnz[i] += rowlengths[j]; 168986d161a7SShri Abhyankar } 169086d161a7SShri Abhyankar } 169186d161a7SShri Abhyankar ierr = PetscFree(rowlengths);CHKERRQ(ierr); 169286d161a7SShri Abhyankar 169386d161a7SShri Abhyankar /* determine max buffer needed and allocate it */ 169486d161a7SShri Abhyankar maxnz = 0; 169586d161a7SShri Abhyankar for (i=0; i<size; i++) { 169686d161a7SShri Abhyankar maxnz = PetscMax(maxnz,procsnz[i]); 169786d161a7SShri Abhyankar } 1698785e854fSJed Brown ierr = PetscMalloc1(maxnz,&cols);CHKERRQ(ierr); 169986d161a7SShri Abhyankar 170086d161a7SShri Abhyankar /* read in my part of the matrix column indices */ 170186d161a7SShri Abhyankar nz = procsnz[0]; 1702785e854fSJed Brown ierr = PetscMalloc1(nz,&mycols);CHKERRQ(ierr); 170386d161a7SShri Abhyankar ierr = PetscBinaryRead(fd,mycols,nz,PETSC_INT);CHKERRQ(ierr); 170486d161a7SShri Abhyankar 170586d161a7SShri Abhyankar /* read in every one elses and ship off */ 170686d161a7SShri Abhyankar for (i=1; i<size; i++) { 170786d161a7SShri Abhyankar nz = procsnz[i]; 170886d161a7SShri Abhyankar ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT);CHKERRQ(ierr); 170986d161a7SShri Abhyankar ierr = MPI_Send(cols,nz,MPIU_INT,i,tag,comm);CHKERRQ(ierr); 171086d161a7SShri Abhyankar } 171186d161a7SShri Abhyankar ierr = PetscFree(cols);CHKERRQ(ierr); 171286d161a7SShri Abhyankar } else { 171386d161a7SShri Abhyankar /* determine buffer space needed for message */ 171486d161a7SShri Abhyankar nz = 0; 171586d161a7SShri Abhyankar for (i=0; i<m; i++) { 171686d161a7SShri Abhyankar nz += ourlens[i]; 171786d161a7SShri Abhyankar } 1718854ce69bSBarry Smith ierr = PetscMalloc1(nz+1,&mycols);CHKERRQ(ierr); 171986d161a7SShri Abhyankar 172086d161a7SShri Abhyankar /* receive message of column indices*/ 172186d161a7SShri Abhyankar ierr = MPI_Recv(mycols,nz,MPIU_INT,0,tag,comm,&status);CHKERRQ(ierr); 172286d161a7SShri Abhyankar ierr = MPI_Get_count(&status,MPIU_INT,&maxnz);CHKERRQ(ierr); 172386d161a7SShri Abhyankar if (maxnz != nz) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"something is wrong with file"); 172486d161a7SShri Abhyankar } 172586d161a7SShri Abhyankar 172649467e7dSBarry Smith ierr = MatSetSizes(newmat,m,n,M,N);CHKERRQ(ierr); 1727dfdea288SBarry Smith a = (Mat_MPIDense*)newmat->data; 17282e3566b0SBarry Smith if (!a->A || !((Mat_SeqDense*)(a->A->data))->user_alloc) { 17290298fd71SBarry Smith ierr = MatMPIDenseSetPreallocation(newmat,NULL);CHKERRQ(ierr); 1730dfdea288SBarry Smith } 173186d161a7SShri Abhyankar 173286d161a7SShri Abhyankar if (!rank) { 1733785e854fSJed Brown ierr = PetscMalloc1(maxnz,&vals);CHKERRQ(ierr); 173486d161a7SShri Abhyankar 173586d161a7SShri Abhyankar /* read in my part of the matrix numerical values */ 173686d161a7SShri Abhyankar nz = procsnz[0]; 173786d161a7SShri Abhyankar ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr); 173886d161a7SShri Abhyankar 173986d161a7SShri Abhyankar /* insert into matrix */ 174086d161a7SShri Abhyankar jj = rstart; 174186d161a7SShri Abhyankar smycols = mycols; 174286d161a7SShri Abhyankar svals = vals; 174386d161a7SShri Abhyankar for (i=0; i<m; i++) { 174486d161a7SShri Abhyankar ierr = MatSetValues(newmat,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr); 174586d161a7SShri Abhyankar smycols += ourlens[i]; 174686d161a7SShri Abhyankar svals += ourlens[i]; 174786d161a7SShri Abhyankar jj++; 174886d161a7SShri Abhyankar } 174986d161a7SShri Abhyankar 175086d161a7SShri Abhyankar /* read in other processors and ship out */ 175186d161a7SShri Abhyankar for (i=1; i<size; i++) { 175286d161a7SShri Abhyankar nz = procsnz[i]; 175386d161a7SShri Abhyankar ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr); 175486d161a7SShri Abhyankar ierr = MPI_Send(vals,nz,MPIU_SCALAR,i,((PetscObject)newmat)->tag,comm);CHKERRQ(ierr); 175586d161a7SShri Abhyankar } 175686d161a7SShri Abhyankar ierr = PetscFree(procsnz);CHKERRQ(ierr); 175786d161a7SShri Abhyankar } else { 175886d161a7SShri Abhyankar /* receive numeric values */ 1759854ce69bSBarry Smith ierr = PetscMalloc1(nz+1,&vals);CHKERRQ(ierr); 176086d161a7SShri Abhyankar 176186d161a7SShri Abhyankar /* receive message of values*/ 176286d161a7SShri Abhyankar ierr = MPI_Recv(vals,nz,MPIU_SCALAR,0,((PetscObject)newmat)->tag,comm,&status);CHKERRQ(ierr); 176386d161a7SShri Abhyankar ierr = MPI_Get_count(&status,MPIU_SCALAR,&maxnz);CHKERRQ(ierr); 176486d161a7SShri Abhyankar if (maxnz != nz) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"something is wrong with file"); 176586d161a7SShri Abhyankar 176686d161a7SShri Abhyankar /* insert into matrix */ 176786d161a7SShri Abhyankar jj = rstart; 176886d161a7SShri Abhyankar smycols = mycols; 176986d161a7SShri Abhyankar svals = vals; 177086d161a7SShri Abhyankar for (i=0; i<m; i++) { 177186d161a7SShri Abhyankar ierr = MatSetValues(newmat,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr); 177286d161a7SShri Abhyankar smycols += ourlens[i]; 177386d161a7SShri Abhyankar svals += ourlens[i]; 177486d161a7SShri Abhyankar jj++; 177586d161a7SShri Abhyankar } 177686d161a7SShri Abhyankar } 177749467e7dSBarry Smith ierr = PetscFree(ourlens);CHKERRQ(ierr); 177886d161a7SShri Abhyankar ierr = PetscFree(vals);CHKERRQ(ierr); 177986d161a7SShri Abhyankar ierr = PetscFree(mycols);CHKERRQ(ierr); 178086d161a7SShri Abhyankar ierr = PetscFree(rowners);CHKERRQ(ierr); 178186d161a7SShri Abhyankar 178286d161a7SShri Abhyankar ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 178386d161a7SShri Abhyankar ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 178486d161a7SShri Abhyankar PetscFunctionReturn(0); 178586d161a7SShri Abhyankar } 178686d161a7SShri Abhyankar 1787ace3abfcSBarry Smith PetscErrorCode MatEqual_MPIDense(Mat A,Mat B,PetscBool *flag) 17886e4ee0c6SHong Zhang { 17896e4ee0c6SHong Zhang Mat_MPIDense *matB = (Mat_MPIDense*)B->data,*matA = (Mat_MPIDense*)A->data; 17906e4ee0c6SHong Zhang Mat a,b; 1791ace3abfcSBarry Smith PetscBool flg; 17926e4ee0c6SHong Zhang PetscErrorCode ierr; 179390ace30eSBarry Smith 17946e4ee0c6SHong Zhang PetscFunctionBegin; 17956e4ee0c6SHong Zhang a = matA->A; 17966e4ee0c6SHong Zhang b = matB->A; 17976e4ee0c6SHong Zhang ierr = MatEqual(a,b,&flg);CHKERRQ(ierr); 1798b2566f29SBarry Smith ierr = MPIU_Allreduce(&flg,flag,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 17996e4ee0c6SHong Zhang PetscFunctionReturn(0); 18006e4ee0c6SHong Zhang } 180190ace30eSBarry Smith 1802baa3c1c6SHong Zhang PetscErrorCode MatDestroy_MatTransMatMult_MPIDense_MPIDense(Mat A) 1803baa3c1c6SHong Zhang { 1804baa3c1c6SHong Zhang PetscErrorCode ierr; 1805baa3c1c6SHong Zhang Mat_MPIDense *a = (Mat_MPIDense*)A->data; 1806baa3c1c6SHong Zhang Mat_TransMatMultDense *atb = a->atbdense; 1807baa3c1c6SHong Zhang 1808baa3c1c6SHong Zhang PetscFunctionBegin; 1809c5ef1628SHong Zhang ierr = PetscFree3(atb->sendbuf,atb->atbarray,atb->recvcounts);CHKERRQ(ierr); 1810baa3c1c6SHong Zhang ierr = (atb->destroy)(A);CHKERRQ(ierr); 1811baa3c1c6SHong Zhang ierr = PetscFree(atb);CHKERRQ(ierr); 1812baa3c1c6SHong Zhang PetscFunctionReturn(0); 1813baa3c1c6SHong Zhang } 1814baa3c1c6SHong Zhang 1815cb20be35SHong Zhang PetscErrorCode MatTransposeMatMultNumeric_MPIDense_MPIDense(Mat A,Mat B,Mat C) 1816cb20be35SHong Zhang { 1817baa3c1c6SHong Zhang Mat_MPIDense *a=(Mat_MPIDense*)A->data, *b=(Mat_MPIDense*)B->data, *c=(Mat_MPIDense*)C->data; 1818660d5466SHong Zhang Mat_SeqDense *aseq=(Mat_SeqDense*)(a->A)->data, *bseq=(Mat_SeqDense*)(b->A)->data; 1819baa3c1c6SHong Zhang Mat_TransMatMultDense *atb = c->atbdense; 1820cb20be35SHong Zhang PetscErrorCode ierr; 1821cb20be35SHong Zhang MPI_Comm comm; 1822d5017740SHong Zhang PetscMPIInt rank,size,*recvcounts=atb->recvcounts; 1823c5ef1628SHong Zhang PetscScalar *carray,*atbarray=atb->atbarray,*sendbuf=atb->sendbuf; 1824d5017740SHong Zhang PetscInt i,cN=C->cmap->N,cM=C->rmap->N,proc,k,j; 1825660d5466SHong Zhang PetscScalar _DOne=1.0,_DZero=0.0; 1826adc7a786SHong Zhang PetscBLASInt am,an,bn,aN; 1827e68c0b26SHong Zhang const PetscInt *ranges; 1828cb20be35SHong Zhang 1829cb20be35SHong Zhang PetscFunctionBegin; 1830cb20be35SHong Zhang ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 1831cb20be35SHong Zhang ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 1832cb20be35SHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1833e68c0b26SHong Zhang 1834c5ef1628SHong Zhang /* compute atbarray = aseq^T * bseq */ 1835660d5466SHong Zhang ierr = PetscBLASIntCast(a->A->cmap->n,&an);CHKERRQ(ierr); 1836660d5466SHong Zhang ierr = PetscBLASIntCast(b->A->cmap->n,&bn);CHKERRQ(ierr); 1837660d5466SHong Zhang ierr = PetscBLASIntCast(a->A->rmap->n,&am);CHKERRQ(ierr); 1838adc7a786SHong Zhang ierr = PetscBLASIntCast(A->cmap->N,&aN);CHKERRQ(ierr); 1839adc7a786SHong Zhang PetscStackCallBLAS("BLASgemm",BLASgemm_("T","N",&an,&bn,&am,&_DOne,aseq->v,&aseq->lda,bseq->v,&bseq->lda,&_DZero,atbarray,&aN)); 1840cb20be35SHong Zhang 1841cb20be35SHong Zhang ierr = MatGetOwnershipRanges(C,&ranges);CHKERRQ(ierr); 1842c5ef1628SHong Zhang for (i=0; i<size; i++) recvcounts[i] = (ranges[i+1] - ranges[i])*cN; 1843cb20be35SHong Zhang 1844660d5466SHong Zhang /* arrange atbarray into sendbuf */ 1845cb20be35SHong Zhang k = 0; 1846cb20be35SHong Zhang for (proc=0; proc<size; proc++) { 1847baa3c1c6SHong Zhang for (j=0; j<cN; j++) { 1848c5ef1628SHong Zhang for (i=ranges[proc]; i<ranges[proc+1]; i++) sendbuf[k++] = atbarray[i+j*cM]; 1849cb20be35SHong Zhang } 1850cb20be35SHong Zhang } 1851c5ef1628SHong Zhang /* sum all atbarray to local values of C */ 1852660d5466SHong Zhang ierr = MatDenseGetArray(c->A,&carray);CHKERRQ(ierr); 18533462b7efSHong Zhang ierr = MPI_Reduce_scatter(sendbuf,carray,recvcounts,MPIU_SCALAR,MPIU_SUM,comm);CHKERRQ(ierr); 1854660d5466SHong Zhang ierr = MatDenseRestoreArray(c->A,&carray);CHKERRQ(ierr); 1855cb20be35SHong Zhang PetscFunctionReturn(0); 1856cb20be35SHong Zhang } 1857cb20be35SHong Zhang 1858cb20be35SHong Zhang PetscErrorCode MatTransposeMatMultSymbolic_MPIDense_MPIDense(Mat A,Mat B,PetscReal fill,Mat *C) 1859cb20be35SHong Zhang { 1860cb20be35SHong Zhang PetscErrorCode ierr; 1861baa3c1c6SHong Zhang Mat Cdense; 1862cb20be35SHong Zhang MPI_Comm comm; 1863baa3c1c6SHong Zhang PetscMPIInt size; 1864660d5466SHong Zhang PetscInt cm=A->cmap->n,cM,cN=B->cmap->N; 1865baa3c1c6SHong Zhang Mat_MPIDense *c; 1866baa3c1c6SHong Zhang Mat_TransMatMultDense *atb; 1867cb20be35SHong Zhang 1868cb20be35SHong Zhang PetscFunctionBegin; 1869baa3c1c6SHong Zhang ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 1870cb20be35SHong Zhang if (A->rmap->rstart != B->rmap->rstart || A->rmap->rend != B->rmap->rend) { 1871cb20be35SHong Zhang SETERRQ4(comm,PETSC_ERR_ARG_SIZ,"Matrix local dimensions are incompatible, A (%D, %D) != B (%D,%D)",A->rmap->rstart,A->rmap->rend,B->rmap->rstart,B->rmap->rend); 1872cb20be35SHong Zhang } 1873cb20be35SHong Zhang 1874baa3c1c6SHong Zhang /* create matrix product Cdense */ 1875baa3c1c6SHong Zhang ierr = MatCreate(comm,&Cdense);CHKERRQ(ierr); 1876660d5466SHong Zhang ierr = MatSetSizes(Cdense,cm,B->cmap->n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr); 1877baa3c1c6SHong Zhang ierr = MatSetType(Cdense,MATMPIDENSE);CHKERRQ(ierr); 1878baa3c1c6SHong Zhang ierr = MatMPIDenseSetPreallocation(Cdense,NULL);CHKERRQ(ierr); 1879baa3c1c6SHong Zhang ierr = MatAssemblyBegin(Cdense,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1880baa3c1c6SHong Zhang ierr = MatAssemblyEnd(Cdense,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1881baa3c1c6SHong Zhang *C = Cdense; 1882baa3c1c6SHong Zhang 1883baa3c1c6SHong Zhang /* create data structure for reuse Cdense */ 1884baa3c1c6SHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1885baa3c1c6SHong Zhang ierr = PetscNew(&atb);CHKERRQ(ierr); 1886660d5466SHong Zhang cM = Cdense->rmap->N; 1887c5ef1628SHong Zhang ierr = PetscMalloc3(cM*cN,&atb->sendbuf,cM*cN,&atb->atbarray,size,&atb->recvcounts);CHKERRQ(ierr); 1888baa3c1c6SHong Zhang 1889baa3c1c6SHong Zhang c = (Mat_MPIDense*)Cdense->data; 1890baa3c1c6SHong Zhang c->atbdense = atb; 1891baa3c1c6SHong Zhang atb->destroy = Cdense->ops->destroy; 1892baa3c1c6SHong Zhang Cdense->ops->destroy = MatDestroy_MatTransMatMult_MPIDense_MPIDense; 1893cb20be35SHong Zhang PetscFunctionReturn(0); 1894cb20be35SHong Zhang } 1895cb20be35SHong Zhang 1896cb20be35SHong Zhang PetscErrorCode MatTransposeMatMult_MPIDense_MPIDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 1897cb20be35SHong Zhang { 1898cb20be35SHong Zhang PetscErrorCode ierr; 1899cb20be35SHong Zhang 1900cb20be35SHong Zhang PetscFunctionBegin; 1901cb20be35SHong Zhang if (scall == MAT_INITIAL_MATRIX) { 1902cb20be35SHong Zhang ierr = MatTransposeMatMultSymbolic_MPIDense_MPIDense(A,B,fill,C);CHKERRQ(ierr); 1903cb20be35SHong Zhang } 1904cb20be35SHong Zhang ierr = MatTransposeMatMultNumeric_MPIDense_MPIDense(A,B,*C);CHKERRQ(ierr); 1905cb20be35SHong Zhang PetscFunctionReturn(0); 1906cb20be35SHong Zhang } 1907320f2790SHong Zhang 1908320f2790SHong Zhang PetscErrorCode MatDestroy_MatMatMult_MPIDense_MPIDense(Mat A) 1909320f2790SHong Zhang { 1910320f2790SHong Zhang PetscErrorCode ierr; 1911320f2790SHong Zhang Mat_MPIDense *a = (Mat_MPIDense*)A->data; 1912320f2790SHong Zhang Mat_MatMultDense *ab = a->abdense; 1913320f2790SHong Zhang 1914320f2790SHong Zhang PetscFunctionBegin; 1915320f2790SHong Zhang ierr = MatDestroy(&ab->Ce);CHKERRQ(ierr); 1916320f2790SHong Zhang ierr = MatDestroy(&ab->Ae);CHKERRQ(ierr); 1917320f2790SHong Zhang ierr = MatDestroy(&ab->Be);CHKERRQ(ierr); 1918320f2790SHong Zhang 1919320f2790SHong Zhang ierr = (ab->destroy)(A);CHKERRQ(ierr); 1920320f2790SHong Zhang ierr = PetscFree(ab);CHKERRQ(ierr); 1921320f2790SHong Zhang PetscFunctionReturn(0); 1922320f2790SHong Zhang } 1923320f2790SHong Zhang 19245242a7b1SHong Zhang #if defined(PETSC_HAVE_ELEMENTAL) 1925320f2790SHong Zhang PetscErrorCode MatMatMultNumeric_MPIDense_MPIDense(Mat A,Mat B,Mat C) 1926320f2790SHong Zhang { 1927320f2790SHong Zhang PetscErrorCode ierr; 1928320f2790SHong Zhang Mat_MPIDense *c=(Mat_MPIDense*)C->data; 1929320f2790SHong Zhang Mat_MatMultDense *ab=c->abdense; 1930320f2790SHong Zhang 1931320f2790SHong Zhang PetscFunctionBegin; 1932de0a22f0SHong Zhang ierr = MatConvert_MPIDense_Elemental(A,MATELEMENTAL,MAT_REUSE_MATRIX, &ab->Ae);CHKERRQ(ierr); 1933de0a22f0SHong Zhang ierr = MatConvert_MPIDense_Elemental(B,MATELEMENTAL,MAT_REUSE_MATRIX, &ab->Be);CHKERRQ(ierr); 1934320f2790SHong Zhang ierr = MatMatMultNumeric(ab->Ae,ab->Be,ab->Ce);CHKERRQ(ierr); 1935de0a22f0SHong Zhang ierr = MatConvert(ab->Ce,MATMPIDENSE,MAT_REUSE_MATRIX,&C);CHKERRQ(ierr); 1936320f2790SHong Zhang PetscFunctionReturn(0); 1937320f2790SHong Zhang } 1938320f2790SHong Zhang 1939320f2790SHong Zhang PetscErrorCode MatMatMultSymbolic_MPIDense_MPIDense(Mat A,Mat B,PetscReal fill,Mat *C) 1940320f2790SHong Zhang { 1941320f2790SHong Zhang PetscErrorCode ierr; 1942320f2790SHong Zhang Mat Ae,Be,Ce; 1943320f2790SHong Zhang Mat_MPIDense *c; 1944320f2790SHong Zhang Mat_MatMultDense *ab; 1945320f2790SHong Zhang 1946320f2790SHong Zhang PetscFunctionBegin; 1947320f2790SHong Zhang if (A->cmap->rstart != B->rmap->rstart || A->cmap->rend != B->rmap->rend) { 1948378336b6SHong Zhang SETERRQ4(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_SIZ,"Matrix local dimensions are incompatible, A (%D, %D) != B (%D,%D)",A->rmap->rstart,A->rmap->rend,B->rmap->rstart,B->rmap->rend); 1949320f2790SHong Zhang } 1950320f2790SHong Zhang 1951320f2790SHong Zhang /* convert A and B to Elemental matrices Ae and Be */ 1952320f2790SHong Zhang ierr = MatConvert(A,MATELEMENTAL,MAT_INITIAL_MATRIX, &Ae);CHKERRQ(ierr); 1953320f2790SHong Zhang ierr = MatConvert(B,MATELEMENTAL,MAT_INITIAL_MATRIX, &Be);CHKERRQ(ierr); 1954320f2790SHong Zhang 1955320f2790SHong Zhang /* Ce = Ae*Be */ 1956320f2790SHong Zhang ierr = MatMatMultSymbolic(Ae,Be,fill,&Ce);CHKERRQ(ierr); 1957320f2790SHong Zhang ierr = MatMatMultNumeric(Ae,Be,Ce);CHKERRQ(ierr); 1958320f2790SHong Zhang 1959320f2790SHong Zhang /* convert Ce to C */ 1960320f2790SHong Zhang ierr = MatConvert(Ce,MATMPIDENSE,MAT_INITIAL_MATRIX,C);CHKERRQ(ierr); 1961320f2790SHong Zhang 1962320f2790SHong Zhang /* create data structure for reuse Cdense */ 1963320f2790SHong Zhang ierr = PetscNew(&ab);CHKERRQ(ierr); 1964320f2790SHong Zhang c = (Mat_MPIDense*)(*C)->data; 1965320f2790SHong Zhang c->abdense = ab; 1966320f2790SHong Zhang 1967320f2790SHong Zhang ab->Ae = Ae; 1968320f2790SHong Zhang ab->Be = Be; 1969320f2790SHong Zhang ab->Ce = Ce; 1970320f2790SHong Zhang ab->destroy = (*C)->ops->destroy; 1971320f2790SHong Zhang (*C)->ops->destroy = MatDestroy_MatMatMult_MPIDense_MPIDense; 19729775878aSHong Zhang (*C)->ops->matmultnumeric = MatMatMultNumeric_MPIDense_MPIDense; 1973320f2790SHong Zhang PetscFunctionReturn(0); 1974320f2790SHong Zhang } 1975320f2790SHong Zhang 1976150d2497SBarry Smith PETSC_INTERN PetscErrorCode MatMatMult_MPIDense_MPIDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 1977320f2790SHong Zhang { 1978320f2790SHong Zhang PetscErrorCode ierr; 1979320f2790SHong Zhang 1980320f2790SHong Zhang PetscFunctionBegin; 1981320f2790SHong Zhang if (scall == MAT_INITIAL_MATRIX) { /* simbolic product includes numeric product */ 198257299f2fSHong Zhang ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr); 1983320f2790SHong Zhang ierr = MatMatMultSymbolic_MPIDense_MPIDense(A,B,fill,C);CHKERRQ(ierr); 198457299f2fSHong Zhang ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr); 1985320f2790SHong Zhang } else { 198657299f2fSHong Zhang ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr); 1987320f2790SHong Zhang ierr = MatMatMultNumeric_MPIDense_MPIDense(A,B,*C);CHKERRQ(ierr); 198857299f2fSHong Zhang ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr); 1989320f2790SHong Zhang } 1990320f2790SHong Zhang PetscFunctionReturn(0); 1991320f2790SHong Zhang } 19925242a7b1SHong Zhang #endif 1993*86aefd0dSHong Zhang 1994