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 147*d3042a70SBarry 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 157*d3042a70SBarry Smith static PetscErrorCode MatDensePlaceArray_MPIDense(Mat A,const PetscScalar array[]) 158*d3042a70SBarry Smith { 159*d3042a70SBarry Smith Mat_MPIDense *a = (Mat_MPIDense*)A->data; 160*d3042a70SBarry Smith PetscErrorCode ierr; 161*d3042a70SBarry Smith 162*d3042a70SBarry Smith PetscFunctionBegin; 163*d3042a70SBarry Smith ierr = MatDensePlaceArray(a->A,array);CHKERRQ(ierr); 164*d3042a70SBarry Smith PetscFunctionReturn(0); 165*d3042a70SBarry Smith } 166*d3042a70SBarry Smith 167*d3042a70SBarry Smith static PetscErrorCode MatDenseResetArray_MPIDense(Mat A) 168*d3042a70SBarry Smith { 169*d3042a70SBarry Smith Mat_MPIDense *a = (Mat_MPIDense*)A->data; 170*d3042a70SBarry Smith PetscErrorCode ierr; 171*d3042a70SBarry Smith 172*d3042a70SBarry Smith PetscFunctionBegin; 173*d3042a70SBarry Smith ierr = MatDenseResetArray(a->A);CHKERRQ(ierr); 174*d3042a70SBarry Smith PetscFunctionReturn(0); 175*d3042a70SBarry Smith } 176*d3042a70SBarry 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); 554*d3042a70SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatDensePlaceArray_C",NULL);CHKERRQ(ierr); 555*d3042a70SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatDenseResetArray_C",NULL);CHKERRQ(ierr); 556*d3042a70SBarry 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); 5683a40ed3dSBarry Smith PetscFunctionReturn(0); 5698965ea79SLois Curfman McInnes } 57039ddd567SLois Curfman McInnes 5716849ba73SBarry Smith static PetscErrorCode MatView_MPIDense_Binary(Mat mat,PetscViewer viewer) 5728965ea79SLois Curfman McInnes { 57339ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 574dfbe8321SBarry Smith PetscErrorCode ierr; 575aa05aa95SBarry Smith PetscViewerFormat format; 576aa05aa95SBarry Smith int fd; 577d0f46423SBarry Smith PetscInt header[4],mmax,N = mat->cmap->N,i,j,m,k; 578aa05aa95SBarry Smith PetscMPIInt rank,tag = ((PetscObject)viewer)->tag,size; 579578230a0SSatish Balay PetscScalar *work,*v,*vv; 580aa05aa95SBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)mdn->A->data; 5817056b6fcSBarry Smith 5823a40ed3dSBarry Smith PetscFunctionBegin; 58339ddd567SLois Curfman McInnes if (mdn->size == 1) { 58439ddd567SLois Curfman McInnes ierr = MatView(mdn->A,viewer);CHKERRQ(ierr); 585aa05aa95SBarry Smith } else { 586aa05aa95SBarry Smith ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 587ce94432eSBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)mat),&rank);CHKERRQ(ierr); 588ce94432eSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)mat),&size);CHKERRQ(ierr); 589aa05aa95SBarry Smith 590aa05aa95SBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 591f4403165SShri Abhyankar if (format == PETSC_VIEWER_NATIVE) { 592aa05aa95SBarry Smith 593aa05aa95SBarry Smith if (!rank) { 594aa05aa95SBarry Smith /* store the matrix as a dense matrix */ 5950700a824SBarry Smith header[0] = MAT_FILE_CLASSID; 596d0f46423SBarry Smith header[1] = mat->rmap->N; 597aa05aa95SBarry Smith header[2] = N; 598aa05aa95SBarry Smith header[3] = MATRIX_BINARY_FORMAT_DENSE; 599aa05aa95SBarry Smith ierr = PetscBinaryWrite(fd,header,4,PETSC_INT,PETSC_TRUE);CHKERRQ(ierr); 600aa05aa95SBarry Smith 601aa05aa95SBarry Smith /* get largest work array needed for transposing array */ 602d0f46423SBarry Smith mmax = mat->rmap->n; 603aa05aa95SBarry Smith for (i=1; i<size; i++) { 604d0f46423SBarry Smith mmax = PetscMax(mmax,mat->rmap->range[i+1] - mat->rmap->range[i]); 6058965ea79SLois Curfman McInnes } 606785e854fSJed Brown ierr = PetscMalloc1(mmax*N,&work);CHKERRQ(ierr); 607aa05aa95SBarry Smith 608aa05aa95SBarry Smith /* write out local array, by rows */ 609d0f46423SBarry Smith m = mat->rmap->n; 610aa05aa95SBarry Smith v = a->v; 611aa05aa95SBarry Smith for (j=0; j<N; j++) { 612aa05aa95SBarry Smith for (i=0; i<m; i++) { 613578230a0SSatish Balay work[j + i*N] = *v++; 614aa05aa95SBarry Smith } 615aa05aa95SBarry Smith } 616aa05aa95SBarry Smith ierr = PetscBinaryWrite(fd,work,m*N,PETSC_SCALAR,PETSC_FALSE);CHKERRQ(ierr); 617aa05aa95SBarry Smith /* get largest work array to receive messages from other processes, excludes process zero */ 618aa05aa95SBarry Smith mmax = 0; 619aa05aa95SBarry Smith for (i=1; i<size; i++) { 620d0f46423SBarry Smith mmax = PetscMax(mmax,mat->rmap->range[i+1] - mat->rmap->range[i]); 621aa05aa95SBarry Smith } 622785e854fSJed Brown ierr = PetscMalloc1(mmax*N,&vv);CHKERRQ(ierr); 623aa05aa95SBarry Smith for (k = 1; k < size; k++) { 624f8009846SMatthew Knepley v = vv; 625d0f46423SBarry Smith m = mat->rmap->range[k+1] - mat->rmap->range[k]; 626ce94432eSBarry Smith ierr = MPIULong_Recv(v,m*N,MPIU_SCALAR,k,tag,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr); 627aa05aa95SBarry Smith 628aa05aa95SBarry Smith for (j = 0; j < N; j++) { 629aa05aa95SBarry Smith for (i = 0; i < m; i++) { 630578230a0SSatish Balay work[j + i*N] = *v++; 631aa05aa95SBarry Smith } 632aa05aa95SBarry Smith } 633aa05aa95SBarry Smith ierr = PetscBinaryWrite(fd,work,m*N,PETSC_SCALAR,PETSC_FALSE);CHKERRQ(ierr); 634aa05aa95SBarry Smith } 635aa05aa95SBarry Smith ierr = PetscFree(work);CHKERRQ(ierr); 636578230a0SSatish Balay ierr = PetscFree(vv);CHKERRQ(ierr); 637aa05aa95SBarry Smith } else { 638ce94432eSBarry Smith ierr = MPIULong_Send(a->v,mat->rmap->n*mat->cmap->N,MPIU_SCALAR,0,tag,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr); 639aa05aa95SBarry Smith } 6406a9046bcSBarry Smith } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"To store a parallel dense matrix you must first call PetscViewerPushFormat(viewer,PETSC_VIEWER_NATIVE)"); 641aa05aa95SBarry Smith } 6423a40ed3dSBarry Smith PetscFunctionReturn(0); 6438965ea79SLois Curfman McInnes } 6448965ea79SLois Curfman McInnes 6457da1fb6eSBarry Smith extern PetscErrorCode MatView_SeqDense(Mat,PetscViewer); 6469804daf3SBarry Smith #include <petscdraw.h> 6476849ba73SBarry Smith static PetscErrorCode MatView_MPIDense_ASCIIorDraworSocket(Mat mat,PetscViewer viewer) 6488965ea79SLois Curfman McInnes { 64939ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 650dfbe8321SBarry Smith PetscErrorCode ierr; 6517da1fb6eSBarry Smith PetscMPIInt rank = mdn->rank; 65219fd82e9SBarry Smith PetscViewerType vtype; 653ace3abfcSBarry Smith PetscBool iascii,isdraw; 654b0a32e0cSBarry Smith PetscViewer sviewer; 655f3ef73ceSBarry Smith PetscViewerFormat format; 6568965ea79SLois Curfman McInnes 6573a40ed3dSBarry Smith PetscFunctionBegin; 658251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 659251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr); 66032077d6dSBarry Smith if (iascii) { 661b0a32e0cSBarry Smith ierr = PetscViewerGetType(viewer,&vtype);CHKERRQ(ierr); 662b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 663456192e2SBarry Smith if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) { 6644e220ebcSLois Curfman McInnes MatInfo info; 665888f2ed8SSatish Balay ierr = MatGetInfo(mat,MAT_LOCAL,&info);CHKERRQ(ierr); 6661575c14dSBarry Smith ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr); 6677b23a99aSBarry 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); 668b0a32e0cSBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 6691575c14dSBarry Smith ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr); 6705d00a290SHong Zhang ierr = VecScatterView(mdn->Mvctx,viewer);CHKERRQ(ierr); 6713a40ed3dSBarry Smith PetscFunctionReturn(0); 672fb9695e5SSatish Balay } else if (format == PETSC_VIEWER_ASCII_INFO) { 6733a40ed3dSBarry Smith PetscFunctionReturn(0); 6748965ea79SLois Curfman McInnes } 675f1af5d2fSBarry Smith } else if (isdraw) { 676b0a32e0cSBarry Smith PetscDraw draw; 677ace3abfcSBarry Smith PetscBool isnull; 678f1af5d2fSBarry Smith 679b0a32e0cSBarry Smith ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 680b0a32e0cSBarry Smith ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); 681f1af5d2fSBarry Smith if (isnull) PetscFunctionReturn(0); 682f1af5d2fSBarry Smith } 68377ed5343SBarry Smith 6847da1fb6eSBarry Smith { 6858965ea79SLois Curfman McInnes /* assemble the entire matrix onto first processor. */ 6868965ea79SLois Curfman McInnes Mat A; 687d0f46423SBarry Smith PetscInt M = mat->rmap->N,N = mat->cmap->N,m,row,i,nz; 688ba8c8a56SBarry Smith PetscInt *cols; 689ba8c8a56SBarry Smith PetscScalar *vals; 6908965ea79SLois Curfman McInnes 691ce94432eSBarry Smith ierr = MatCreate(PetscObjectComm((PetscObject)mat),&A);CHKERRQ(ierr); 6928965ea79SLois Curfman McInnes if (!rank) { 693f69a0ea3SMatthew Knepley ierr = MatSetSizes(A,M,N,M,N);CHKERRQ(ierr); 6943a40ed3dSBarry Smith } else { 695f69a0ea3SMatthew Knepley ierr = MatSetSizes(A,0,0,M,N);CHKERRQ(ierr); 6968965ea79SLois Curfman McInnes } 6977adad957SLisandro Dalcin /* Since this is a temporary matrix, MATMPIDENSE instead of ((PetscObject)A)->type_name here is probably acceptable. */ 698878740d9SKris Buschelman ierr = MatSetType(A,MATMPIDENSE);CHKERRQ(ierr); 6990298fd71SBarry Smith ierr = MatMPIDenseSetPreallocation(A,NULL);CHKERRQ(ierr); 7003bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)A);CHKERRQ(ierr); 7018965ea79SLois Curfman McInnes 70239ddd567SLois Curfman McInnes /* Copy the matrix ... This isn't the most efficient means, 70339ddd567SLois Curfman McInnes but it's quick for now */ 70451022da4SBarry Smith A->insertmode = INSERT_VALUES; 7052205254eSKarl Rupp 7062205254eSKarl Rupp row = mat->rmap->rstart; 7072205254eSKarl Rupp m = mdn->A->rmap->n; 7088965ea79SLois Curfman McInnes for (i=0; i<m; i++) { 709ba8c8a56SBarry Smith ierr = MatGetRow_MPIDense(mat,row,&nz,&cols,&vals);CHKERRQ(ierr); 710ba8c8a56SBarry Smith ierr = MatSetValues_MPIDense(A,1,&row,nz,cols,vals,INSERT_VALUES);CHKERRQ(ierr); 711ba8c8a56SBarry Smith ierr = MatRestoreRow_MPIDense(mat,row,&nz,&cols,&vals);CHKERRQ(ierr); 71239ddd567SLois Curfman McInnes row++; 7138965ea79SLois Curfman McInnes } 7148965ea79SLois Curfman McInnes 7156d4a8577SBarry Smith ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 7166d4a8577SBarry Smith ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 7173f08860eSBarry Smith ierr = PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr); 718b9b97703SBarry Smith if (!rank) { 7191a9d3c3cSBarry Smith ierr = PetscObjectSetName((PetscObject)((Mat_MPIDense*)(A->data))->A,((PetscObject)mat)->name);CHKERRQ(ierr); 7207da1fb6eSBarry Smith ierr = MatView_SeqDense(((Mat_MPIDense*)(A->data))->A,sviewer);CHKERRQ(ierr); 7218965ea79SLois Curfman McInnes } 7223f08860eSBarry Smith ierr = PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr); 723b0a32e0cSBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 7246bf464f9SBarry Smith ierr = MatDestroy(&A);CHKERRQ(ierr); 7258965ea79SLois Curfman McInnes } 7263a40ed3dSBarry Smith PetscFunctionReturn(0); 7278965ea79SLois Curfman McInnes } 7288965ea79SLois Curfman McInnes 729dfbe8321SBarry Smith PetscErrorCode MatView_MPIDense(Mat mat,PetscViewer viewer) 7308965ea79SLois Curfman McInnes { 731dfbe8321SBarry Smith PetscErrorCode ierr; 732ace3abfcSBarry Smith PetscBool iascii,isbinary,isdraw,issocket; 7338965ea79SLois Curfman McInnes 734433994e6SBarry Smith PetscFunctionBegin; 735251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 736251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr); 737251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSOCKET,&issocket);CHKERRQ(ierr); 738251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr); 7390f5bd95cSBarry Smith 74032077d6dSBarry Smith if (iascii || issocket || isdraw) { 741f1af5d2fSBarry Smith ierr = MatView_MPIDense_ASCIIorDraworSocket(mat,viewer);CHKERRQ(ierr); 7420f5bd95cSBarry Smith } else if (isbinary) { 7433a40ed3dSBarry Smith ierr = MatView_MPIDense_Binary(mat,viewer);CHKERRQ(ierr); 74411aeaf0aSBarry Smith } 7453a40ed3dSBarry Smith PetscFunctionReturn(0); 7468965ea79SLois Curfman McInnes } 7478965ea79SLois Curfman McInnes 748dfbe8321SBarry Smith PetscErrorCode MatGetInfo_MPIDense(Mat A,MatInfoType flag,MatInfo *info) 7498965ea79SLois Curfman McInnes { 7503501a2bdSLois Curfman McInnes Mat_MPIDense *mat = (Mat_MPIDense*)A->data; 7513501a2bdSLois Curfman McInnes Mat mdn = mat->A; 752dfbe8321SBarry Smith PetscErrorCode ierr; 753329f5518SBarry Smith PetscReal isend[5],irecv[5]; 7548965ea79SLois Curfman McInnes 7553a40ed3dSBarry Smith PetscFunctionBegin; 7564e220ebcSLois Curfman McInnes info->block_size = 1.0; 7572205254eSKarl Rupp 7584e220ebcSLois Curfman McInnes ierr = MatGetInfo(mdn,MAT_LOCAL,info);CHKERRQ(ierr); 7592205254eSKarl Rupp 7604e220ebcSLois Curfman McInnes isend[0] = info->nz_used; isend[1] = info->nz_allocated; isend[2] = info->nz_unneeded; 7614e220ebcSLois Curfman McInnes isend[3] = info->memory; isend[4] = info->mallocs; 7628965ea79SLois Curfman McInnes if (flag == MAT_LOCAL) { 7634e220ebcSLois Curfman McInnes info->nz_used = isend[0]; 7644e220ebcSLois Curfman McInnes info->nz_allocated = isend[1]; 7654e220ebcSLois Curfman McInnes info->nz_unneeded = isend[2]; 7664e220ebcSLois Curfman McInnes info->memory = isend[3]; 7674e220ebcSLois Curfman McInnes info->mallocs = isend[4]; 7688965ea79SLois Curfman McInnes } else if (flag == MAT_GLOBAL_MAX) { 769b2566f29SBarry Smith ierr = MPIU_Allreduce(isend,irecv,5,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 7702205254eSKarl Rupp 7714e220ebcSLois Curfman McInnes info->nz_used = irecv[0]; 7724e220ebcSLois Curfman McInnes info->nz_allocated = irecv[1]; 7734e220ebcSLois Curfman McInnes info->nz_unneeded = irecv[2]; 7744e220ebcSLois Curfman McInnes info->memory = irecv[3]; 7754e220ebcSLois Curfman McInnes info->mallocs = irecv[4]; 7768965ea79SLois Curfman McInnes } else if (flag == MAT_GLOBAL_SUM) { 777b2566f29SBarry Smith ierr = MPIU_Allreduce(isend,irecv,5,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 7782205254eSKarl Rupp 7794e220ebcSLois Curfman McInnes info->nz_used = irecv[0]; 7804e220ebcSLois Curfman McInnes info->nz_allocated = irecv[1]; 7814e220ebcSLois Curfman McInnes info->nz_unneeded = irecv[2]; 7824e220ebcSLois Curfman McInnes info->memory = irecv[3]; 7834e220ebcSLois Curfman McInnes info->mallocs = irecv[4]; 7848965ea79SLois Curfman McInnes } 7854e220ebcSLois Curfman McInnes info->fill_ratio_given = 0; /* no parallel LU/ILU/Cholesky */ 7864e220ebcSLois Curfman McInnes info->fill_ratio_needed = 0; 7874e220ebcSLois Curfman McInnes info->factor_mallocs = 0; 7883a40ed3dSBarry Smith PetscFunctionReturn(0); 7898965ea79SLois Curfman McInnes } 7908965ea79SLois Curfman McInnes 791ace3abfcSBarry Smith PetscErrorCode MatSetOption_MPIDense(Mat A,MatOption op,PetscBool flg) 7928965ea79SLois Curfman McInnes { 79339ddd567SLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense*)A->data; 794dfbe8321SBarry Smith PetscErrorCode ierr; 7958965ea79SLois Curfman McInnes 7963a40ed3dSBarry Smith PetscFunctionBegin; 79712c028f9SKris Buschelman switch (op) { 798512a5fc5SBarry Smith case MAT_NEW_NONZERO_LOCATIONS: 79912c028f9SKris Buschelman case MAT_NEW_NONZERO_LOCATION_ERR: 80012c028f9SKris Buschelman case MAT_NEW_NONZERO_ALLOCATION_ERR: 80143674050SBarry Smith MatCheckPreallocated(A,1); 8024e0d8c25SBarry Smith ierr = MatSetOption(a->A,op,flg);CHKERRQ(ierr); 80312c028f9SKris Buschelman break; 80412c028f9SKris Buschelman case MAT_ROW_ORIENTED: 80543674050SBarry Smith MatCheckPreallocated(A,1); 8064e0d8c25SBarry Smith a->roworiented = flg; 8074e0d8c25SBarry Smith ierr = MatSetOption(a->A,op,flg);CHKERRQ(ierr); 80812c028f9SKris Buschelman break; 8094e0d8c25SBarry Smith case MAT_NEW_DIAGONALS: 81013fa8e87SLisandro Dalcin case MAT_KEEP_NONZERO_PATTERN: 81112c028f9SKris Buschelman case MAT_USE_HASH_TABLE: 812290bbb0aSBarry Smith ierr = PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);CHKERRQ(ierr); 81312c028f9SKris Buschelman break; 81412c028f9SKris Buschelman case MAT_IGNORE_OFF_PROC_ENTRIES: 8154e0d8c25SBarry Smith a->donotstash = flg; 81612c028f9SKris Buschelman break; 81777e54ba9SKris Buschelman case MAT_SYMMETRIC: 81877e54ba9SKris Buschelman case MAT_STRUCTURALLY_SYMMETRIC: 8199a4540c5SBarry Smith case MAT_HERMITIAN: 8209a4540c5SBarry Smith case MAT_SYMMETRY_ETERNAL: 821600fe468SBarry Smith case MAT_IGNORE_LOWER_TRIANGULAR: 822290bbb0aSBarry Smith ierr = PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);CHKERRQ(ierr); 82377e54ba9SKris Buschelman break; 82412c028f9SKris Buschelman default: 825e32f2f54SBarry Smith SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"unknown option %s",MatOptions[op]); 8263a40ed3dSBarry Smith } 8273a40ed3dSBarry Smith PetscFunctionReturn(0); 8288965ea79SLois Curfman McInnes } 8298965ea79SLois Curfman McInnes 8308965ea79SLois Curfman McInnes 831dfbe8321SBarry Smith PetscErrorCode MatDiagonalScale_MPIDense(Mat A,Vec ll,Vec rr) 8325b2fa520SLois Curfman McInnes { 8335b2fa520SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)A->data; 8345b2fa520SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)mdn->A->data; 83587828ca2SBarry Smith PetscScalar *l,*r,x,*v; 836dfbe8321SBarry Smith PetscErrorCode ierr; 837d0f46423SBarry Smith PetscInt i,j,s2a,s3a,s2,s3,m=mdn->A->rmap->n,n=mdn->A->cmap->n; 8385b2fa520SLois Curfman McInnes 8395b2fa520SLois Curfman McInnes PetscFunctionBegin; 84072d926a5SLois Curfman McInnes ierr = MatGetLocalSize(A,&s2,&s3);CHKERRQ(ierr); 8415b2fa520SLois Curfman McInnes if (ll) { 84272d926a5SLois Curfman McInnes ierr = VecGetLocalSize(ll,&s2a);CHKERRQ(ierr); 843e32f2f54SBarry Smith if (s2a != s2) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Left scaling vector non-conforming local size, %d != %d.", s2a, s2); 8441ebc52fbSHong Zhang ierr = VecGetArray(ll,&l);CHKERRQ(ierr); 8455b2fa520SLois Curfman McInnes for (i=0; i<m; i++) { 8465b2fa520SLois Curfman McInnes x = l[i]; 8475b2fa520SLois Curfman McInnes v = mat->v + i; 8485b2fa520SLois Curfman McInnes for (j=0; j<n; j++) { (*v) *= x; v+= m;} 8495b2fa520SLois Curfman McInnes } 8501ebc52fbSHong Zhang ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr); 851efee365bSSatish Balay ierr = PetscLogFlops(n*m);CHKERRQ(ierr); 8525b2fa520SLois Curfman McInnes } 8535b2fa520SLois Curfman McInnes if (rr) { 854175be7b4SMatthew Knepley ierr = VecGetLocalSize(rr,&s3a);CHKERRQ(ierr); 855e32f2f54SBarry Smith if (s3a != s3) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Right scaling vec non-conforming local size, %d != %d.", s3a, s3); 856ca9f406cSSatish Balay ierr = VecScatterBegin(mdn->Mvctx,rr,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 857ca9f406cSSatish Balay ierr = VecScatterEnd(mdn->Mvctx,rr,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 8581ebc52fbSHong Zhang ierr = VecGetArray(mdn->lvec,&r);CHKERRQ(ierr); 8595b2fa520SLois Curfman McInnes for (i=0; i<n; i++) { 8605b2fa520SLois Curfman McInnes x = r[i]; 8615b2fa520SLois Curfman McInnes v = mat->v + i*m; 8622205254eSKarl Rupp for (j=0; j<m; j++) (*v++) *= x; 8635b2fa520SLois Curfman McInnes } 8641ebc52fbSHong Zhang ierr = VecRestoreArray(mdn->lvec,&r);CHKERRQ(ierr); 865efee365bSSatish Balay ierr = PetscLogFlops(n*m);CHKERRQ(ierr); 8665b2fa520SLois Curfman McInnes } 8675b2fa520SLois Curfman McInnes PetscFunctionReturn(0); 8685b2fa520SLois Curfman McInnes } 8695b2fa520SLois Curfman McInnes 870dfbe8321SBarry Smith PetscErrorCode MatNorm_MPIDense(Mat A,NormType type,PetscReal *nrm) 871096963f5SLois Curfman McInnes { 8723501a2bdSLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)A->data; 8733501a2bdSLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)mdn->A->data; 874dfbe8321SBarry Smith PetscErrorCode ierr; 87513f74950SBarry Smith PetscInt i,j; 876329f5518SBarry Smith PetscReal sum = 0.0; 87787828ca2SBarry Smith PetscScalar *v = mat->v; 8783501a2bdSLois Curfman McInnes 8793a40ed3dSBarry Smith PetscFunctionBegin; 8803501a2bdSLois Curfman McInnes if (mdn->size == 1) { 881064f8208SBarry Smith ierr = MatNorm(mdn->A,type,nrm);CHKERRQ(ierr); 8823501a2bdSLois Curfman McInnes } else { 8833501a2bdSLois Curfman McInnes if (type == NORM_FROBENIUS) { 884d0f46423SBarry Smith for (i=0; i<mdn->A->cmap->n*mdn->A->rmap->n; i++) { 885329f5518SBarry Smith sum += PetscRealPart(PetscConj(*v)*(*v)); v++; 8863501a2bdSLois Curfman McInnes } 887b2566f29SBarry Smith ierr = MPIU_Allreduce(&sum,nrm,1,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 8888f1a2a5eSBarry Smith *nrm = PetscSqrtReal(*nrm); 889dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*mdn->A->cmap->n*mdn->A->rmap->n);CHKERRQ(ierr); 8903a40ed3dSBarry Smith } else if (type == NORM_1) { 891329f5518SBarry Smith PetscReal *tmp,*tmp2; 892dcca6d9dSJed Brown ierr = PetscMalloc2(A->cmap->N,&tmp,A->cmap->N,&tmp2);CHKERRQ(ierr); 89374ed9c26SBarry Smith ierr = PetscMemzero(tmp,A->cmap->N*sizeof(PetscReal));CHKERRQ(ierr); 89474ed9c26SBarry Smith ierr = PetscMemzero(tmp2,A->cmap->N*sizeof(PetscReal));CHKERRQ(ierr); 895064f8208SBarry Smith *nrm = 0.0; 8963501a2bdSLois Curfman McInnes v = mat->v; 897d0f46423SBarry Smith for (j=0; j<mdn->A->cmap->n; j++) { 898d0f46423SBarry Smith for (i=0; i<mdn->A->rmap->n; i++) { 89967e560aaSBarry Smith tmp[j] += PetscAbsScalar(*v); v++; 9003501a2bdSLois Curfman McInnes } 9013501a2bdSLois Curfman McInnes } 902b2566f29SBarry Smith ierr = MPIU_Allreduce(tmp,tmp2,A->cmap->N,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 903d0f46423SBarry Smith for (j=0; j<A->cmap->N; j++) { 904064f8208SBarry Smith if (tmp2[j] > *nrm) *nrm = tmp2[j]; 9053501a2bdSLois Curfman McInnes } 9068627564fSBarry Smith ierr = PetscFree2(tmp,tmp2);CHKERRQ(ierr); 907d0f46423SBarry Smith ierr = PetscLogFlops(A->cmap->n*A->rmap->n);CHKERRQ(ierr); 9083a40ed3dSBarry Smith } else if (type == NORM_INFINITY) { /* max row norm */ 909329f5518SBarry Smith PetscReal ntemp; 9103501a2bdSLois Curfman McInnes ierr = MatNorm(mdn->A,type,&ntemp);CHKERRQ(ierr); 911b2566f29SBarry Smith ierr = MPIU_Allreduce(&ntemp,nrm,1,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 912ce94432eSBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support for two norm"); 9133501a2bdSLois Curfman McInnes } 9143a40ed3dSBarry Smith PetscFunctionReturn(0); 9153501a2bdSLois Curfman McInnes } 9163501a2bdSLois Curfman McInnes 917fc4dec0aSBarry Smith PetscErrorCode MatTranspose_MPIDense(Mat A,MatReuse reuse,Mat *matout) 9183501a2bdSLois Curfman McInnes { 9193501a2bdSLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense*)A->data; 9203501a2bdSLois Curfman McInnes Mat_SeqDense *Aloc = (Mat_SeqDense*)a->A->data; 9213501a2bdSLois Curfman McInnes Mat B; 922d0f46423SBarry Smith PetscInt M = A->rmap->N,N = A->cmap->N,m,n,*rwork,rstart = A->rmap->rstart; 9236849ba73SBarry Smith PetscErrorCode ierr; 92413f74950SBarry Smith PetscInt j,i; 92587828ca2SBarry Smith PetscScalar *v; 9263501a2bdSLois Curfman McInnes 9273a40ed3dSBarry Smith PetscFunctionBegin; 928cf37664fSBarry Smith if (reuse == MAT_INPLACE_MATRIX && M != N) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Supports square matrix only in-place"); 929cf37664fSBarry Smith if (reuse == MAT_INITIAL_MATRIX || reuse == MAT_INPLACE_MATRIX) { 930ce94432eSBarry Smith ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr); 931d0f46423SBarry Smith ierr = MatSetSizes(B,A->cmap->n,A->rmap->n,N,M);CHKERRQ(ierr); 9327adad957SLisandro Dalcin ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr); 9330298fd71SBarry Smith ierr = MatMPIDenseSetPreallocation(B,NULL);CHKERRQ(ierr); 934fc4dec0aSBarry Smith } else { 935fc4dec0aSBarry Smith B = *matout; 936fc4dec0aSBarry Smith } 9373501a2bdSLois Curfman McInnes 938d0f46423SBarry Smith m = a->A->rmap->n; n = a->A->cmap->n; v = Aloc->v; 939785e854fSJed Brown ierr = PetscMalloc1(m,&rwork);CHKERRQ(ierr); 9403501a2bdSLois Curfman McInnes for (i=0; i<m; i++) rwork[i] = rstart + i; 9411acff37aSSatish Balay for (j=0; j<n; j++) { 9423501a2bdSLois Curfman McInnes ierr = MatSetValues(B,1,&j,m,rwork,v,INSERT_VALUES);CHKERRQ(ierr); 9433501a2bdSLois Curfman McInnes v += m; 9443501a2bdSLois Curfman McInnes } 945606d414cSSatish Balay ierr = PetscFree(rwork);CHKERRQ(ierr); 9466d4a8577SBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 9476d4a8577SBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 948cf37664fSBarry Smith if (reuse == MAT_INITIAL_MATRIX || reuse == MAT_REUSE_MATRIX) { 9493501a2bdSLois Curfman McInnes *matout = B; 9503501a2bdSLois Curfman McInnes } else { 95128be2f97SBarry Smith ierr = MatHeaderMerge(A,&B);CHKERRQ(ierr); 9523501a2bdSLois Curfman McInnes } 9533a40ed3dSBarry Smith PetscFunctionReturn(0); 954096963f5SLois Curfman McInnes } 955096963f5SLois Curfman McInnes 95644cd7ae7SLois Curfman McInnes 9576849ba73SBarry Smith static PetscErrorCode MatDuplicate_MPIDense(Mat,MatDuplicateOption,Mat*); 958d84338a6SBarry Smith extern PetscErrorCode MatScale_MPIDense(Mat,PetscScalar); 9598965ea79SLois Curfman McInnes 9604994cf47SJed Brown PetscErrorCode MatSetUp_MPIDense(Mat A) 961273d9f13SBarry Smith { 962dfbe8321SBarry Smith PetscErrorCode ierr; 963273d9f13SBarry Smith 964273d9f13SBarry Smith PetscFunctionBegin; 965273d9f13SBarry Smith ierr = MatMPIDenseSetPreallocation(A,0);CHKERRQ(ierr); 966273d9f13SBarry Smith PetscFunctionReturn(0); 967273d9f13SBarry Smith } 968273d9f13SBarry Smith 969488007eeSBarry Smith PetscErrorCode MatAXPY_MPIDense(Mat Y,PetscScalar alpha,Mat X,MatStructure str) 970488007eeSBarry Smith { 971488007eeSBarry Smith PetscErrorCode ierr; 972488007eeSBarry Smith Mat_MPIDense *A = (Mat_MPIDense*)Y->data, *B = (Mat_MPIDense*)X->data; 973488007eeSBarry Smith 974488007eeSBarry Smith PetscFunctionBegin; 975488007eeSBarry Smith ierr = MatAXPY(A->A,alpha,B->A,str);CHKERRQ(ierr); 976a3fa217bSJose E. Roman ierr = PetscObjectStateIncrease((PetscObject)Y);CHKERRQ(ierr); 977488007eeSBarry Smith PetscFunctionReturn(0); 978488007eeSBarry Smith } 979488007eeSBarry Smith 9807087cfbeSBarry Smith PetscErrorCode MatConjugate_MPIDense(Mat mat) 981ba337c44SJed Brown { 982ba337c44SJed Brown Mat_MPIDense *a = (Mat_MPIDense*)mat->data; 983ba337c44SJed Brown PetscErrorCode ierr; 984ba337c44SJed Brown 985ba337c44SJed Brown PetscFunctionBegin; 986ba337c44SJed Brown ierr = MatConjugate(a->A);CHKERRQ(ierr); 987ba337c44SJed Brown PetscFunctionReturn(0); 988ba337c44SJed Brown } 989ba337c44SJed Brown 990ba337c44SJed Brown PetscErrorCode MatRealPart_MPIDense(Mat A) 991ba337c44SJed Brown { 992ba337c44SJed Brown Mat_MPIDense *a = (Mat_MPIDense*)A->data; 993ba337c44SJed Brown PetscErrorCode ierr; 994ba337c44SJed Brown 995ba337c44SJed Brown PetscFunctionBegin; 996ba337c44SJed Brown ierr = MatRealPart(a->A);CHKERRQ(ierr); 997ba337c44SJed Brown PetscFunctionReturn(0); 998ba337c44SJed Brown } 999ba337c44SJed Brown 1000ba337c44SJed Brown PetscErrorCode MatImaginaryPart_MPIDense(Mat A) 1001ba337c44SJed Brown { 1002ba337c44SJed Brown Mat_MPIDense *a = (Mat_MPIDense*)A->data; 1003ba337c44SJed Brown PetscErrorCode ierr; 1004ba337c44SJed Brown 1005ba337c44SJed Brown PetscFunctionBegin; 1006ba337c44SJed Brown ierr = MatImaginaryPart(a->A);CHKERRQ(ierr); 1007ba337c44SJed Brown PetscFunctionReturn(0); 1008ba337c44SJed Brown } 1009ba337c44SJed Brown 10100716a85fSBarry Smith extern PetscErrorCode MatGetColumnNorms_SeqDense(Mat,NormType,PetscReal*); 10110716a85fSBarry Smith PetscErrorCode MatGetColumnNorms_MPIDense(Mat A,NormType type,PetscReal *norms) 10120716a85fSBarry Smith { 10130716a85fSBarry Smith PetscErrorCode ierr; 10140716a85fSBarry Smith PetscInt i,n; 10150716a85fSBarry Smith Mat_MPIDense *a = (Mat_MPIDense*) A->data; 10160716a85fSBarry Smith PetscReal *work; 10170716a85fSBarry Smith 10180716a85fSBarry Smith PetscFunctionBegin; 10190298fd71SBarry Smith ierr = MatGetSize(A,NULL,&n);CHKERRQ(ierr); 1020785e854fSJed Brown ierr = PetscMalloc1(n,&work);CHKERRQ(ierr); 10210716a85fSBarry Smith ierr = MatGetColumnNorms_SeqDense(a->A,type,work);CHKERRQ(ierr); 10220716a85fSBarry Smith if (type == NORM_2) { 10230716a85fSBarry Smith for (i=0; i<n; i++) work[i] *= work[i]; 10240716a85fSBarry Smith } 10250716a85fSBarry Smith if (type == NORM_INFINITY) { 1026b2566f29SBarry Smith ierr = MPIU_Allreduce(work,norms,n,MPIU_REAL,MPIU_MAX,A->hdr.comm);CHKERRQ(ierr); 10270716a85fSBarry Smith } else { 1028b2566f29SBarry Smith ierr = MPIU_Allreduce(work,norms,n,MPIU_REAL,MPIU_SUM,A->hdr.comm);CHKERRQ(ierr); 10290716a85fSBarry Smith } 10300716a85fSBarry Smith ierr = PetscFree(work);CHKERRQ(ierr); 10310716a85fSBarry Smith if (type == NORM_2) { 10328f1a2a5eSBarry Smith for (i=0; i<n; i++) norms[i] = PetscSqrtReal(norms[i]); 10330716a85fSBarry Smith } 10340716a85fSBarry Smith PetscFunctionReturn(0); 10350716a85fSBarry Smith } 10360716a85fSBarry Smith 103773a71a0fSBarry Smith static PetscErrorCode MatSetRandom_MPIDense(Mat x,PetscRandom rctx) 103873a71a0fSBarry Smith { 103973a71a0fSBarry Smith Mat_MPIDense *d = (Mat_MPIDense*)x->data; 104073a71a0fSBarry Smith PetscErrorCode ierr; 104173a71a0fSBarry Smith PetscScalar *a; 104273a71a0fSBarry Smith PetscInt m,n,i; 104373a71a0fSBarry Smith 104473a71a0fSBarry Smith PetscFunctionBegin; 104573a71a0fSBarry Smith ierr = MatGetSize(d->A,&m,&n);CHKERRQ(ierr); 10468c778c55SBarry Smith ierr = MatDenseGetArray(d->A,&a);CHKERRQ(ierr); 104773a71a0fSBarry Smith for (i=0; i<m*n; i++) { 104873a71a0fSBarry Smith ierr = PetscRandomGetValue(rctx,a+i);CHKERRQ(ierr); 104973a71a0fSBarry Smith } 10508c778c55SBarry Smith ierr = MatDenseRestoreArray(d->A,&a);CHKERRQ(ierr); 105173a71a0fSBarry Smith PetscFunctionReturn(0); 105273a71a0fSBarry Smith } 105373a71a0fSBarry Smith 1054fd4e9aacSBarry Smith extern PetscErrorCode MatMatMultNumeric_MPIDense(Mat A,Mat,Mat); 1055fd4e9aacSBarry Smith 10563b49f96aSBarry Smith static PetscErrorCode MatMissingDiagonal_MPIDense(Mat A,PetscBool *missing,PetscInt *d) 10573b49f96aSBarry Smith { 10583b49f96aSBarry Smith PetscFunctionBegin; 10593b49f96aSBarry Smith *missing = PETSC_FALSE; 10603b49f96aSBarry Smith PetscFunctionReturn(0); 10613b49f96aSBarry Smith } 10623b49f96aSBarry Smith 10638965ea79SLois Curfman McInnes /* -------------------------------------------------------------------*/ 106409dc0095SBarry Smith static struct _MatOps MatOps_Values = { MatSetValues_MPIDense, 106509dc0095SBarry Smith MatGetRow_MPIDense, 106609dc0095SBarry Smith MatRestoreRow_MPIDense, 106709dc0095SBarry Smith MatMult_MPIDense, 106897304618SKris Buschelman /* 4*/ MatMultAdd_MPIDense, 10697c922b88SBarry Smith MatMultTranspose_MPIDense, 10707c922b88SBarry Smith MatMultTransposeAdd_MPIDense, 10718965ea79SLois Curfman McInnes 0, 107209dc0095SBarry Smith 0, 107309dc0095SBarry Smith 0, 107497304618SKris Buschelman /* 10*/ 0, 107509dc0095SBarry Smith 0, 107609dc0095SBarry Smith 0, 107709dc0095SBarry Smith 0, 107809dc0095SBarry Smith MatTranspose_MPIDense, 107997304618SKris Buschelman /* 15*/ MatGetInfo_MPIDense, 10806e4ee0c6SHong Zhang MatEqual_MPIDense, 108109dc0095SBarry Smith MatGetDiagonal_MPIDense, 10825b2fa520SLois Curfman McInnes MatDiagonalScale_MPIDense, 108309dc0095SBarry Smith MatNorm_MPIDense, 108497304618SKris Buschelman /* 20*/ MatAssemblyBegin_MPIDense, 108509dc0095SBarry Smith MatAssemblyEnd_MPIDense, 108609dc0095SBarry Smith MatSetOption_MPIDense, 108709dc0095SBarry Smith MatZeroEntries_MPIDense, 1088d519adbfSMatthew Knepley /* 24*/ MatZeroRows_MPIDense, 1089919b68f7SBarry Smith 0, 109001b82886SBarry Smith 0, 109101b82886SBarry Smith 0, 109201b82886SBarry Smith 0, 10934994cf47SJed Brown /* 29*/ MatSetUp_MPIDense, 1094273d9f13SBarry Smith 0, 109509dc0095SBarry Smith 0, 1096c56a70eeSBarry Smith MatGetDiagonalBlock_MPIDense, 10978c778c55SBarry Smith 0, 1098d519adbfSMatthew Knepley /* 34*/ MatDuplicate_MPIDense, 109909dc0095SBarry Smith 0, 110009dc0095SBarry Smith 0, 110109dc0095SBarry Smith 0, 110209dc0095SBarry Smith 0, 1103d519adbfSMatthew Knepley /* 39*/ MatAXPY_MPIDense, 11047dae84e0SHong Zhang MatCreateSubMatrices_MPIDense, 110509dc0095SBarry Smith 0, 110609dc0095SBarry Smith MatGetValues_MPIDense, 110709dc0095SBarry Smith 0, 1108d519adbfSMatthew Knepley /* 44*/ 0, 110909dc0095SBarry Smith MatScale_MPIDense, 11107d68702bSBarry Smith MatShift_Basic, 111109dc0095SBarry Smith 0, 111209dc0095SBarry Smith 0, 111373a71a0fSBarry Smith /* 49*/ MatSetRandom_MPIDense, 111409dc0095SBarry Smith 0, 111509dc0095SBarry Smith 0, 111609dc0095SBarry Smith 0, 111709dc0095SBarry Smith 0, 1118d519adbfSMatthew Knepley /* 54*/ 0, 111909dc0095SBarry Smith 0, 112009dc0095SBarry Smith 0, 112109dc0095SBarry Smith 0, 112209dc0095SBarry Smith 0, 11237dae84e0SHong Zhang /* 59*/ MatCreateSubMatrix_MPIDense, 1124b9b97703SBarry Smith MatDestroy_MPIDense, 1125b9b97703SBarry Smith MatView_MPIDense, 1126357abbc8SBarry Smith 0, 112797304618SKris Buschelman 0, 1128d519adbfSMatthew Knepley /* 64*/ 0, 112997304618SKris Buschelman 0, 113097304618SKris Buschelman 0, 113197304618SKris Buschelman 0, 113297304618SKris Buschelman 0, 1133d519adbfSMatthew Knepley /* 69*/ 0, 113497304618SKris Buschelman 0, 113597304618SKris Buschelman 0, 113697304618SKris Buschelman 0, 113797304618SKris Buschelman 0, 1138d519adbfSMatthew Knepley /* 74*/ 0, 113997304618SKris Buschelman 0, 114097304618SKris Buschelman 0, 114197304618SKris Buschelman 0, 114297304618SKris Buschelman 0, 1143d519adbfSMatthew Knepley /* 79*/ 0, 114497304618SKris Buschelman 0, 114597304618SKris Buschelman 0, 114697304618SKris Buschelman 0, 11475bba2384SShri Abhyankar /* 83*/ MatLoad_MPIDense, 1148865e5f61SKris Buschelman 0, 1149865e5f61SKris Buschelman 0, 1150865e5f61SKris Buschelman 0, 1151865e5f61SKris Buschelman 0, 1152865e5f61SKris Buschelman 0, 11539775878aSHong Zhang #if defined(PETSC_HAVE_ELEMENTAL) 1154320f2790SHong Zhang /* 89*/ MatMatMult_MPIDense_MPIDense, 1155320f2790SHong Zhang MatMatMultSymbolic_MPIDense_MPIDense, 11569775878aSHong Zhang #else 11579775878aSHong Zhang /* 89*/ 0, 1158865e5f61SKris Buschelman 0, 11599775878aSHong Zhang #endif 1160fd4e9aacSBarry Smith MatMatMultNumeric_MPIDense, 11612fbe02b9SBarry Smith 0, 1162ba337c44SJed Brown 0, 1163d519adbfSMatthew Knepley /* 94*/ 0, 1164865e5f61SKris Buschelman 0, 1165865e5f61SKris Buschelman 0, 1166ba337c44SJed Brown 0, 1167ba337c44SJed Brown 0, 1168ba337c44SJed Brown /* 99*/ 0, 1169ba337c44SJed Brown 0, 1170ba337c44SJed Brown 0, 1171ba337c44SJed Brown MatConjugate_MPIDense, 1172ba337c44SJed Brown 0, 1173ba337c44SJed Brown /*104*/ 0, 1174ba337c44SJed Brown MatRealPart_MPIDense, 1175ba337c44SJed Brown MatImaginaryPart_MPIDense, 117686d161a7SShri Abhyankar 0, 117786d161a7SShri Abhyankar 0, 117886d161a7SShri Abhyankar /*109*/ 0, 117986d161a7SShri Abhyankar 0, 118086d161a7SShri Abhyankar 0, 118186d161a7SShri Abhyankar 0, 11823b49f96aSBarry Smith MatMissingDiagonal_MPIDense, 118386d161a7SShri Abhyankar /*114*/ 0, 118486d161a7SShri Abhyankar 0, 118586d161a7SShri Abhyankar 0, 118686d161a7SShri Abhyankar 0, 118786d161a7SShri Abhyankar 0, 118886d161a7SShri Abhyankar /*119*/ 0, 118986d161a7SShri Abhyankar 0, 119086d161a7SShri Abhyankar 0, 11910716a85fSBarry Smith 0, 11920716a85fSBarry Smith 0, 11930716a85fSBarry Smith /*124*/ 0, 11943964eb88SJed Brown MatGetColumnNorms_MPIDense, 11953964eb88SJed Brown 0, 11963964eb88SJed Brown 0, 11973964eb88SJed Brown 0, 11983964eb88SJed Brown /*129*/ 0, 1199cb20be35SHong Zhang MatTransposeMatMult_MPIDense_MPIDense, 1200cb20be35SHong Zhang MatTransposeMatMultSymbolic_MPIDense_MPIDense, 1201cb20be35SHong Zhang MatTransposeMatMultNumeric_MPIDense_MPIDense, 12023964eb88SJed Brown 0, 12033964eb88SJed Brown /*134*/ 0, 12043964eb88SJed Brown 0, 12053964eb88SJed Brown 0, 12063964eb88SJed Brown 0, 12073964eb88SJed Brown 0, 12083964eb88SJed Brown /*139*/ 0, 1209f9426fe0SMark Adams 0, 12103964eb88SJed Brown 0 1211ba337c44SJed Brown }; 12128965ea79SLois Curfman McInnes 12137087cfbeSBarry Smith PetscErrorCode MatMPIDenseSetPreallocation_MPIDense(Mat mat,PetscScalar *data) 1214a23d5eceSKris Buschelman { 1215a23d5eceSKris Buschelman Mat_MPIDense *a; 1216dfbe8321SBarry Smith PetscErrorCode ierr; 1217a23d5eceSKris Buschelman 1218a23d5eceSKris Buschelman PetscFunctionBegin; 1219a23d5eceSKris Buschelman mat->preallocated = PETSC_TRUE; 1220a23d5eceSKris Buschelman /* Note: For now, when data is specified above, this assumes the user correctly 1221a23d5eceSKris Buschelman allocates the local dense storage space. We should add error checking. */ 1222a23d5eceSKris Buschelman 1223a23d5eceSKris Buschelman a = (Mat_MPIDense*)mat->data; 122434ef9618SShri Abhyankar ierr = PetscLayoutSetUp(mat->rmap);CHKERRQ(ierr); 122534ef9618SShri Abhyankar ierr = PetscLayoutSetUp(mat->cmap);CHKERRQ(ierr); 122634ef9618SShri Abhyankar a->nvec = mat->cmap->n; 122734ef9618SShri Abhyankar 1228f69a0ea3SMatthew Knepley ierr = MatCreate(PETSC_COMM_SELF,&a->A);CHKERRQ(ierr); 1229d0f46423SBarry Smith ierr = MatSetSizes(a->A,mat->rmap->n,mat->cmap->N,mat->rmap->n,mat->cmap->N);CHKERRQ(ierr); 12305c5985e7SKris Buschelman ierr = MatSetType(a->A,MATSEQDENSE);CHKERRQ(ierr); 12315c5985e7SKris Buschelman ierr = MatSeqDenseSetPreallocation(a->A,data);CHKERRQ(ierr); 12323bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)a->A);CHKERRQ(ierr); 1233a23d5eceSKris Buschelman PetscFunctionReturn(0); 1234a23d5eceSKris Buschelman } 1235a23d5eceSKris Buschelman 123665b80a83SHong Zhang #if defined(PETSC_HAVE_ELEMENTAL) 1237cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatConvert_MPIDense_Elemental(Mat A, MatType newtype,MatReuse reuse,Mat *newmat) 12388baccfbdSHong Zhang { 12398ea901baSHong Zhang Mat mat_elemental; 12408ea901baSHong Zhang PetscErrorCode ierr; 124132d7a744SHong Zhang PetscScalar *v; 124232d7a744SHong Zhang PetscInt m=A->rmap->n,N=A->cmap->N,rstart=A->rmap->rstart,i,*rows,*cols; 12438ea901baSHong Zhang 12448baccfbdSHong Zhang PetscFunctionBegin; 1245378336b6SHong Zhang if (reuse == MAT_REUSE_MATRIX) { 1246378336b6SHong Zhang mat_elemental = *newmat; 1247378336b6SHong Zhang ierr = MatZeroEntries(*newmat);CHKERRQ(ierr); 1248378336b6SHong Zhang } else { 1249378336b6SHong Zhang ierr = MatCreate(PetscObjectComm((PetscObject)A), &mat_elemental);CHKERRQ(ierr); 1250378336b6SHong Zhang ierr = MatSetSizes(mat_elemental,PETSC_DECIDE,PETSC_DECIDE,A->rmap->N,A->cmap->N);CHKERRQ(ierr); 1251378336b6SHong Zhang ierr = MatSetType(mat_elemental,MATELEMENTAL);CHKERRQ(ierr); 1252378336b6SHong Zhang ierr = MatSetUp(mat_elemental);CHKERRQ(ierr); 125332d7a744SHong Zhang ierr = MatSetOption(mat_elemental,MAT_ROW_ORIENTED,PETSC_FALSE);CHKERRQ(ierr); 1254378336b6SHong Zhang } 1255378336b6SHong Zhang 125632d7a744SHong Zhang ierr = PetscMalloc2(m,&rows,N,&cols);CHKERRQ(ierr); 125732d7a744SHong Zhang for (i=0; i<N; i++) cols[i] = i; 125832d7a744SHong Zhang for (i=0; i<m; i++) rows[i] = rstart + i; 12598ea901baSHong Zhang 12608ea901baSHong Zhang /* PETSc-Elemental interaface uses axpy for setting off-processor entries, only ADD_VALUES is allowed */ 126132d7a744SHong Zhang ierr = MatDenseGetArray(A,&v);CHKERRQ(ierr); 126232d7a744SHong Zhang ierr = MatSetValues(mat_elemental,m,rows,N,cols,v,ADD_VALUES);CHKERRQ(ierr); 12638ea901baSHong Zhang ierr = MatAssemblyBegin(mat_elemental, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 12648ea901baSHong Zhang ierr = MatAssemblyEnd(mat_elemental, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 126532d7a744SHong Zhang ierr = MatDenseRestoreArray(A,&v);CHKERRQ(ierr); 126632d7a744SHong Zhang ierr = PetscFree2(rows,cols);CHKERRQ(ierr); 12678ea901baSHong Zhang 1268511c6705SHong Zhang if (reuse == MAT_INPLACE_MATRIX) { 126928be2f97SBarry Smith ierr = MatHeaderReplace(A,&mat_elemental);CHKERRQ(ierr); 12708ea901baSHong Zhang } else { 12718ea901baSHong Zhang *newmat = mat_elemental; 12728ea901baSHong Zhang } 12738baccfbdSHong Zhang PetscFunctionReturn(0); 12748baccfbdSHong Zhang } 127565b80a83SHong Zhang #endif 12768baccfbdSHong Zhang 12778cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_MPIDense(Mat mat) 1278273d9f13SBarry Smith { 1279273d9f13SBarry Smith Mat_MPIDense *a; 1280dfbe8321SBarry Smith PetscErrorCode ierr; 1281273d9f13SBarry Smith 1282273d9f13SBarry Smith PetscFunctionBegin; 1283b00a9115SJed Brown ierr = PetscNewLog(mat,&a);CHKERRQ(ierr); 1284b0a32e0cSBarry Smith mat->data = (void*)a; 1285273d9f13SBarry Smith ierr = PetscMemcpy(mat->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 1286273d9f13SBarry Smith 1287273d9f13SBarry Smith mat->insertmode = NOT_SET_VALUES; 1288ce94432eSBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)mat),&a->rank);CHKERRQ(ierr); 1289ce94432eSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)mat),&a->size);CHKERRQ(ierr); 1290273d9f13SBarry Smith 1291273d9f13SBarry Smith /* build cache for off array entries formed */ 1292273d9f13SBarry Smith a->donotstash = PETSC_FALSE; 12932205254eSKarl Rupp 1294ce94432eSBarry Smith ierr = MatStashCreate_Private(PetscObjectComm((PetscObject)mat),1,&mat->stash);CHKERRQ(ierr); 1295273d9f13SBarry Smith 1296273d9f13SBarry Smith /* stuff used for matrix vector multiply */ 1297273d9f13SBarry Smith a->lvec = 0; 1298273d9f13SBarry Smith a->Mvctx = 0; 1299273d9f13SBarry Smith a->roworiented = PETSC_TRUE; 1300273d9f13SBarry Smith 1301bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatDenseGetArray_C",MatDenseGetArray_MPIDense);CHKERRQ(ierr); 1302*d3042a70SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatDensePlaceArray_C",MatDensePlaceArray_MPIDense);CHKERRQ(ierr); 1303*d3042a70SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatDenseResetArray_C",MatDenseResetArray_MPIDense);CHKERRQ(ierr); 1304bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatDenseRestoreArray_C",MatDenseRestoreArray_MPIDense);CHKERRQ(ierr); 13058baccfbdSHong Zhang #if defined(PETSC_HAVE_ELEMENTAL) 13068baccfbdSHong Zhang ierr = PetscObjectComposeFunction((PetscObject)mat,"MatConvert_mpidense_elemental_C",MatConvert_MPIDense_Elemental);CHKERRQ(ierr); 13078baccfbdSHong Zhang #endif 1308bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMPIDenseSetPreallocation_C",MatMPIDenseSetPreallocation_MPIDense);CHKERRQ(ierr); 1309bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMatMult_mpiaij_mpidense_C",MatMatMult_MPIAIJ_MPIDense);CHKERRQ(ierr); 1310bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMatMultSymbolic_mpiaij_mpidense_C",MatMatMultSymbolic_MPIAIJ_MPIDense);CHKERRQ(ierr); 1311bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMatMultNumeric_mpiaij_mpidense_C",MatMatMultNumeric_MPIAIJ_MPIDense);CHKERRQ(ierr); 13128949adfdSHong Zhang 13138949adfdSHong Zhang ierr = PetscObjectComposeFunction((PetscObject)mat,"MatTransposeMatMult_mpiaij_mpidense_C",MatTransposeMatMult_MPIAIJ_MPIDense);CHKERRQ(ierr); 13148949adfdSHong Zhang ierr = PetscObjectComposeFunction((PetscObject)mat,"MatTransposeMatMultSymbolic_mpiaij_mpidense_C",MatTransposeMatMultSymbolic_MPIAIJ_MPIDense);CHKERRQ(ierr); 13158949adfdSHong Zhang ierr = PetscObjectComposeFunction((PetscObject)mat,"MatTransposeMatMultNumeric_mpiaij_mpidense_C",MatTransposeMatMultNumeric_MPIAIJ_MPIDense);CHKERRQ(ierr); 131638aed534SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)mat,MATMPIDENSE);CHKERRQ(ierr); 1317273d9f13SBarry Smith PetscFunctionReturn(0); 1318273d9f13SBarry Smith } 1319273d9f13SBarry Smith 1320209238afSKris Buschelman /*MC 1321002d173eSKris Buschelman MATDENSE - MATDENSE = "dense" - A matrix type to be used for dense matrices. 1322209238afSKris Buschelman 1323209238afSKris Buschelman This matrix type is identical to MATSEQDENSE when constructed with a single process communicator, 1324209238afSKris Buschelman and MATMPIDENSE otherwise. 1325209238afSKris Buschelman 1326209238afSKris Buschelman Options Database Keys: 1327209238afSKris Buschelman . -mat_type dense - sets the matrix type to "dense" during a call to MatSetFromOptions() 1328209238afSKris Buschelman 1329209238afSKris Buschelman Level: beginner 1330209238afSKris Buschelman 133101b82886SBarry Smith 1332209238afSKris Buschelman .seealso: MatCreateMPIDense,MATSEQDENSE,MATMPIDENSE 1333209238afSKris Buschelman M*/ 1334209238afSKris Buschelman 1335273d9f13SBarry Smith /*@C 1336273d9f13SBarry Smith MatMPIDenseSetPreallocation - Sets the array used to store the matrix entries 1337273d9f13SBarry Smith 1338273d9f13SBarry Smith Not collective 1339273d9f13SBarry Smith 1340273d9f13SBarry Smith Input Parameters: 13411c4f3114SJed Brown . B - the matrix 13420298fd71SBarry Smith - data - optional location of matrix data. Set data=NULL for PETSc 1343273d9f13SBarry Smith to control all matrix memory allocation. 1344273d9f13SBarry Smith 1345273d9f13SBarry Smith Notes: 1346273d9f13SBarry Smith The dense format is fully compatible with standard Fortran 77 1347273d9f13SBarry Smith storage by columns. 1348273d9f13SBarry Smith 1349273d9f13SBarry Smith The data input variable is intended primarily for Fortran programmers 1350273d9f13SBarry Smith who wish to allocate their own matrix memory space. Most users should 13510298fd71SBarry Smith set data=NULL. 1352273d9f13SBarry Smith 1353273d9f13SBarry Smith Level: intermediate 1354273d9f13SBarry Smith 1355273d9f13SBarry Smith .keywords: matrix,dense, parallel 1356273d9f13SBarry Smith 1357273d9f13SBarry Smith .seealso: MatCreate(), MatCreateSeqDense(), MatSetValues() 1358273d9f13SBarry Smith @*/ 13591c4f3114SJed Brown PetscErrorCode MatMPIDenseSetPreallocation(Mat B,PetscScalar *data) 1360273d9f13SBarry Smith { 13614ac538c5SBarry Smith PetscErrorCode ierr; 1362273d9f13SBarry Smith 1363273d9f13SBarry Smith PetscFunctionBegin; 13641c4f3114SJed Brown ierr = PetscTryMethod(B,"MatMPIDenseSetPreallocation_C",(Mat,PetscScalar*),(B,data));CHKERRQ(ierr); 1365273d9f13SBarry Smith PetscFunctionReturn(0); 1366273d9f13SBarry Smith } 1367273d9f13SBarry Smith 1368*d3042a70SBarry Smith /*@ 1369*d3042a70SBarry Smith MatDensePlaceArray - Allows one to replace the array in a dense array with an 1370*d3042a70SBarry Smith array provided by the user. This is useful to avoid copying an array 1371*d3042a70SBarry Smith into a matrix 1372*d3042a70SBarry Smith 1373*d3042a70SBarry Smith Not Collective 1374*d3042a70SBarry Smith 1375*d3042a70SBarry Smith Input Parameters: 1376*d3042a70SBarry Smith + mat - the matrix 1377*d3042a70SBarry Smith - array - the array in column major order 1378*d3042a70SBarry Smith 1379*d3042a70SBarry Smith Notes: 1380*d3042a70SBarry 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 1381*d3042a70SBarry Smith freed when the matrix is destroyed. 1382*d3042a70SBarry Smith 1383*d3042a70SBarry Smith Level: developer 1384*d3042a70SBarry Smith 1385*d3042a70SBarry Smith .seealso: MatDenseGetArray(), MatDenseResetArray(), VecPlaceArray(), VecGetArray(), VecRestoreArray(), VecReplaceArray(), VecResetArray() 1386*d3042a70SBarry Smith 1387*d3042a70SBarry Smith @*/ 1388*d3042a70SBarry Smith PetscErrorCode MatDensePlaceArray(Mat mat,const PetscScalar array[]) 1389*d3042a70SBarry Smith { 1390*d3042a70SBarry Smith PetscErrorCode ierr; 1391*d3042a70SBarry Smith PetscFunctionBegin; 1392*d3042a70SBarry Smith ierr = PetscUseMethod(mat,"MatDensePlaceArray_C",(Mat,const PetscScalar*),(mat,array));CHKERRQ(ierr); 1393*d3042a70SBarry Smith ierr = PetscObjectStateIncrease((PetscObject)mat);CHKERRQ(ierr); 1394*d3042a70SBarry Smith PetscFunctionReturn(0); 1395*d3042a70SBarry Smith } 1396*d3042a70SBarry Smith 1397*d3042a70SBarry Smith /*@ 1398*d3042a70SBarry Smith MatDenseResetArray - Resets the matrix array to that it previously had before the call to MatDensePlaceArray() 1399*d3042a70SBarry Smith 1400*d3042a70SBarry Smith Not Collective 1401*d3042a70SBarry Smith 1402*d3042a70SBarry Smith Input Parameters: 1403*d3042a70SBarry Smith . mat - the matrix 1404*d3042a70SBarry Smith 1405*d3042a70SBarry Smith Notes: 1406*d3042a70SBarry Smith You can only call this after a call to MatDensePlaceArray() 1407*d3042a70SBarry Smith 1408*d3042a70SBarry Smith Level: developer 1409*d3042a70SBarry Smith 1410*d3042a70SBarry Smith .seealso: MatDenseGetArray(), MatDensePlaceArray(), VecPlaceArray(), VecGetArray(), VecRestoreArray(), VecReplaceArray(), VecResetArray() 1411*d3042a70SBarry Smith 1412*d3042a70SBarry Smith @*/ 1413*d3042a70SBarry Smith PetscErrorCode MatDenseResetArray(Mat mat) 1414*d3042a70SBarry Smith { 1415*d3042a70SBarry Smith PetscErrorCode ierr; 1416*d3042a70SBarry Smith PetscFunctionBegin; 1417*d3042a70SBarry Smith ierr = PetscUseMethod(mat,"MatDenseResetArray_C",(Mat),(mat));CHKERRQ(ierr); 1418*d3042a70SBarry Smith ierr = PetscObjectStateIncrease((PetscObject)mat);CHKERRQ(ierr); 1419*d3042a70SBarry Smith PetscFunctionReturn(0); 1420*d3042a70SBarry Smith } 1421*d3042a70SBarry Smith 14228965ea79SLois Curfman McInnes /*@C 142369b1f4b7SBarry Smith MatCreateDense - Creates a parallel matrix in dense format. 14248965ea79SLois Curfman McInnes 1425db81eaa0SLois Curfman McInnes Collective on MPI_Comm 1426db81eaa0SLois Curfman McInnes 14278965ea79SLois Curfman McInnes Input Parameters: 1428db81eaa0SLois Curfman McInnes + comm - MPI communicator 14298965ea79SLois Curfman McInnes . m - number of local rows (or PETSC_DECIDE to have calculated if M is given) 1430db81eaa0SLois Curfman McInnes . n - number of local columns (or PETSC_DECIDE to have calculated if N is given) 14318965ea79SLois Curfman McInnes . M - number of global rows (or PETSC_DECIDE to have calculated if m is given) 1432db81eaa0SLois Curfman McInnes . N - number of global columns (or PETSC_DECIDE to have calculated if n is given) 14336cfe35ddSJose E. Roman - data - optional location of matrix data. Set data=NULL (PETSC_NULL_SCALAR for Fortran users) for PETSc 1434dfc5480cSLois Curfman McInnes to control all matrix memory allocation. 14358965ea79SLois Curfman McInnes 14368965ea79SLois Curfman McInnes Output Parameter: 1437477f1c0bSLois Curfman McInnes . A - the matrix 14388965ea79SLois Curfman McInnes 1439b259b22eSLois Curfman McInnes Notes: 144039ddd567SLois Curfman McInnes The dense format is fully compatible with standard Fortran 77 144139ddd567SLois Curfman McInnes storage by columns. 14428965ea79SLois Curfman McInnes 144318f449edSLois Curfman McInnes The data input variable is intended primarily for Fortran programmers 144418f449edSLois Curfman McInnes who wish to allocate their own matrix memory space. Most users should 14456cfe35ddSJose E. Roman set data=NULL (PETSC_NULL_SCALAR for Fortran users). 144618f449edSLois Curfman McInnes 14478965ea79SLois Curfman McInnes The user MUST specify either the local or global matrix dimensions 14488965ea79SLois Curfman McInnes (possibly both). 14498965ea79SLois Curfman McInnes 1450027ccd11SLois Curfman McInnes Level: intermediate 1451027ccd11SLois Curfman McInnes 145239ddd567SLois Curfman McInnes .keywords: matrix,dense, parallel 14538965ea79SLois Curfman McInnes 145439ddd567SLois Curfman McInnes .seealso: MatCreate(), MatCreateSeqDense(), MatSetValues() 14558965ea79SLois Curfman McInnes @*/ 145669b1f4b7SBarry Smith PetscErrorCode MatCreateDense(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscScalar *data,Mat *A) 14578965ea79SLois Curfman McInnes { 14586849ba73SBarry Smith PetscErrorCode ierr; 145913f74950SBarry Smith PetscMPIInt size; 14608965ea79SLois Curfman McInnes 14613a40ed3dSBarry Smith PetscFunctionBegin; 1462f69a0ea3SMatthew Knepley ierr = MatCreate(comm,A);CHKERRQ(ierr); 1463f69a0ea3SMatthew Knepley ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr); 1464273d9f13SBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1465273d9f13SBarry Smith if (size > 1) { 1466273d9f13SBarry Smith ierr = MatSetType(*A,MATMPIDENSE);CHKERRQ(ierr); 1467273d9f13SBarry Smith ierr = MatMPIDenseSetPreallocation(*A,data);CHKERRQ(ierr); 14686cfe35ddSJose E. Roman if (data) { /* user provided data array, so no need to assemble */ 14696cfe35ddSJose E. Roman ierr = MatSetUpMultiply_MPIDense(*A);CHKERRQ(ierr); 14706cfe35ddSJose E. Roman (*A)->assembled = PETSC_TRUE; 14716cfe35ddSJose E. Roman } 1472273d9f13SBarry Smith } else { 1473273d9f13SBarry Smith ierr = MatSetType(*A,MATSEQDENSE);CHKERRQ(ierr); 1474273d9f13SBarry Smith ierr = MatSeqDenseSetPreallocation(*A,data);CHKERRQ(ierr); 14758c469469SLois Curfman McInnes } 14763a40ed3dSBarry Smith PetscFunctionReturn(0); 14778965ea79SLois Curfman McInnes } 14788965ea79SLois Curfman McInnes 14796849ba73SBarry Smith static PetscErrorCode MatDuplicate_MPIDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat) 14808965ea79SLois Curfman McInnes { 14818965ea79SLois Curfman McInnes Mat mat; 14823501a2bdSLois Curfman McInnes Mat_MPIDense *a,*oldmat = (Mat_MPIDense*)A->data; 1483dfbe8321SBarry Smith PetscErrorCode ierr; 14848965ea79SLois Curfman McInnes 14853a40ed3dSBarry Smith PetscFunctionBegin; 14868965ea79SLois Curfman McInnes *newmat = 0; 1487ce94432eSBarry Smith ierr = MatCreate(PetscObjectComm((PetscObject)A),&mat);CHKERRQ(ierr); 1488d0f46423SBarry Smith ierr = MatSetSizes(mat,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr); 14897adad957SLisandro Dalcin ierr = MatSetType(mat,((PetscObject)A)->type_name);CHKERRQ(ierr); 1490834f8fabSBarry Smith a = (Mat_MPIDense*)mat->data; 1491e04c1aa4SHong Zhang ierr = PetscMemcpy(mat->ops,A->ops,sizeof(struct _MatOps));CHKERRQ(ierr); 14925aa7edbeSHong Zhang 1493d5f3da31SBarry Smith mat->factortype = A->factortype; 1494c456f294SBarry Smith mat->assembled = PETSC_TRUE; 1495273d9f13SBarry Smith mat->preallocated = PETSC_TRUE; 14968965ea79SLois Curfman McInnes 14978965ea79SLois Curfman McInnes a->size = oldmat->size; 14988965ea79SLois Curfman McInnes a->rank = oldmat->rank; 1499e0fa3b82SLois Curfman McInnes mat->insertmode = NOT_SET_VALUES; 1500b9b97703SBarry Smith a->nvec = oldmat->nvec; 15013782ba37SSatish Balay a->donotstash = oldmat->donotstash; 1502e04c1aa4SHong Zhang 15031e1e43feSBarry Smith ierr = PetscLayoutReference(A->rmap,&mat->rmap);CHKERRQ(ierr); 15041e1e43feSBarry Smith ierr = PetscLayoutReference(A->cmap,&mat->cmap);CHKERRQ(ierr); 15058965ea79SLois Curfman McInnes 1506329f5518SBarry Smith ierr = MatSetUpMultiply_MPIDense(mat);CHKERRQ(ierr); 15075609ef8eSBarry Smith ierr = MatDuplicate(oldmat->A,cpvalues,&a->A);CHKERRQ(ierr); 15083bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)a->A);CHKERRQ(ierr); 150901b82886SBarry Smith 15108965ea79SLois Curfman McInnes *newmat = mat; 15113a40ed3dSBarry Smith PetscFunctionReturn(0); 15128965ea79SLois Curfman McInnes } 15138965ea79SLois Curfman McInnes 15149d36ed5fSBarry Smith PetscErrorCode MatLoad_MPIDense_DenseInFile(MPI_Comm comm,PetscInt fd,PetscInt M,PetscInt N,Mat newmat) 151586d161a7SShri Abhyankar { 151686d161a7SShri Abhyankar PetscErrorCode ierr; 151786d161a7SShri Abhyankar PetscMPIInt rank,size; 151874c13d95SBarry Smith const PetscInt *rowners; 151974c13d95SBarry Smith PetscInt i,m,n,nz,j,mMax; 152086d161a7SShri Abhyankar PetscScalar *array,*vals,*vals_ptr; 15219d36ed5fSBarry Smith Mat_MPIDense *a = (Mat_MPIDense*)newmat->data; 152286d161a7SShri Abhyankar 152386d161a7SShri Abhyankar PetscFunctionBegin; 152486d161a7SShri Abhyankar ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 152586d161a7SShri Abhyankar ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 152686d161a7SShri Abhyankar 152774c13d95SBarry Smith /* determine ownership of rows and columns */ 152874c13d95SBarry Smith m = (newmat->rmap->n < 0) ? PETSC_DECIDE : newmat->rmap->n; 152974c13d95SBarry Smith n = (newmat->cmap->n < 0) ? PETSC_DECIDE : newmat->cmap->n; 153086d161a7SShri Abhyankar 153174c13d95SBarry Smith ierr = MatSetSizes(newmat,m,n,M,N);CHKERRQ(ierr); 15329d36ed5fSBarry Smith if (!a->A || !((Mat_SeqDense*)(a->A->data))->user_alloc) { 15330298fd71SBarry Smith ierr = MatMPIDenseSetPreallocation(newmat,NULL);CHKERRQ(ierr); 15349d36ed5fSBarry Smith } 15358c778c55SBarry Smith ierr = MatDenseGetArray(newmat,&array);CHKERRQ(ierr); 153649467e7dSBarry Smith ierr = MatGetLocalSize(newmat,&m,NULL);CHKERRQ(ierr); 153774c13d95SBarry Smith ierr = MatGetOwnershipRanges(newmat,&rowners);CHKERRQ(ierr); 1538396e826eSBarry Smith ierr = MPI_Reduce(&m,&mMax,1,MPIU_INT,MPI_MAX,0,comm);CHKERRQ(ierr); 153986d161a7SShri Abhyankar if (!rank) { 15409d36ed5fSBarry Smith ierr = PetscMalloc1(mMax*N,&vals);CHKERRQ(ierr); 154186d161a7SShri Abhyankar 154286d161a7SShri Abhyankar /* read in my part of the matrix numerical values */ 154386d161a7SShri Abhyankar ierr = PetscBinaryRead(fd,vals,m*N,PETSC_SCALAR);CHKERRQ(ierr); 154486d161a7SShri Abhyankar 154586d161a7SShri Abhyankar /* insert into matrix-by row (this is why cannot directly read into array */ 154686d161a7SShri Abhyankar vals_ptr = vals; 154786d161a7SShri Abhyankar for (i=0; i<m; i++) { 154886d161a7SShri Abhyankar for (j=0; j<N; j++) { 154986d161a7SShri Abhyankar array[i + j*m] = *vals_ptr++; 155086d161a7SShri Abhyankar } 155186d161a7SShri Abhyankar } 155286d161a7SShri Abhyankar 155386d161a7SShri Abhyankar /* read in other processors and ship out */ 155486d161a7SShri Abhyankar for (i=1; i<size; i++) { 155586d161a7SShri Abhyankar nz = (rowners[i+1] - rowners[i])*N; 155686d161a7SShri Abhyankar ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr); 1557a25532f0SBarry Smith ierr = MPIULong_Send(vals,nz,MPIU_SCALAR,i,((PetscObject)(newmat))->tag,comm);CHKERRQ(ierr); 155886d161a7SShri Abhyankar } 155986d161a7SShri Abhyankar } else { 156086d161a7SShri Abhyankar /* receive numeric values */ 1561785e854fSJed Brown ierr = PetscMalloc1(m*N,&vals);CHKERRQ(ierr); 156286d161a7SShri Abhyankar 156386d161a7SShri Abhyankar /* receive message of values*/ 1564a25532f0SBarry Smith ierr = MPIULong_Recv(vals,m*N,MPIU_SCALAR,0,((PetscObject)(newmat))->tag,comm);CHKERRQ(ierr); 156586d161a7SShri Abhyankar 156686d161a7SShri Abhyankar /* insert into matrix-by row (this is why cannot directly read into array */ 156786d161a7SShri Abhyankar vals_ptr = vals; 156886d161a7SShri Abhyankar for (i=0; i<m; i++) { 156986d161a7SShri Abhyankar for (j=0; j<N; j++) { 157086d161a7SShri Abhyankar array[i + j*m] = *vals_ptr++; 157186d161a7SShri Abhyankar } 157286d161a7SShri Abhyankar } 157386d161a7SShri Abhyankar } 15748c778c55SBarry Smith ierr = MatDenseRestoreArray(newmat,&array);CHKERRQ(ierr); 157586d161a7SShri Abhyankar ierr = PetscFree(vals);CHKERRQ(ierr); 157686d161a7SShri Abhyankar ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 157786d161a7SShri Abhyankar ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 157886d161a7SShri Abhyankar PetscFunctionReturn(0); 157986d161a7SShri Abhyankar } 158086d161a7SShri Abhyankar 1581112444f4SShri Abhyankar PetscErrorCode MatLoad_MPIDense(Mat newmat,PetscViewer viewer) 158286d161a7SShri Abhyankar { 1583dfdea288SBarry Smith Mat_MPIDense *a; 158486d161a7SShri Abhyankar PetscScalar *vals,*svals; 1585ce94432eSBarry Smith MPI_Comm comm; 158686d161a7SShri Abhyankar MPI_Status status; 158749467e7dSBarry Smith PetscMPIInt rank,size,tag = ((PetscObject)viewer)->tag,*rowners,*sndcounts,m,n,maxnz; 158886d161a7SShri Abhyankar PetscInt header[4],*rowlengths = 0,M,N,*cols; 158949467e7dSBarry Smith PetscInt *ourlens,*procsnz = 0,jj,*mycols,*smycols; 15909d36ed5fSBarry Smith PetscInt i,nz,j,rstart,rend; 159186d161a7SShri Abhyankar int fd; 159286d161a7SShri Abhyankar PetscErrorCode ierr; 159386d161a7SShri Abhyankar 159486d161a7SShri Abhyankar PetscFunctionBegin; 1595c98fd787SBarry Smith /* force binary viewer to load .info file if it has not yet done so */ 1596c98fd787SBarry Smith ierr = PetscViewerSetUp(viewer);CHKERRQ(ierr); 1597ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)viewer,&comm);CHKERRQ(ierr); 159886d161a7SShri Abhyankar ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 159986d161a7SShri Abhyankar ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 160086d161a7SShri Abhyankar ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 16015872f025SBarry Smith if (!rank) { 160286d161a7SShri Abhyankar ierr = PetscBinaryRead(fd,(char*)header,4,PETSC_INT);CHKERRQ(ierr); 160386d161a7SShri Abhyankar if (header[0] != MAT_FILE_CLASSID) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"not matrix object"); 160486d161a7SShri Abhyankar } 160586d161a7SShri Abhyankar ierr = MPI_Bcast(header+1,3,MPIU_INT,0,comm);CHKERRQ(ierr); 160686d161a7SShri Abhyankar M = header[1]; N = header[2]; nz = header[3]; 160786d161a7SShri Abhyankar 160886d161a7SShri Abhyankar /* If global rows/cols are set to PETSC_DECIDE, set it to the sizes given in the file */ 16099d36ed5fSBarry Smith if (newmat->rmap->N < 0) newmat->rmap->N = M; 16109d36ed5fSBarry Smith if (newmat->cmap->N < 0) newmat->cmap->N = N; 161186d161a7SShri Abhyankar 16129d36ed5fSBarry 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); 16139d36ed5fSBarry 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); 161486d161a7SShri Abhyankar 161586d161a7SShri Abhyankar /* 161686d161a7SShri Abhyankar Handle case where matrix is stored on disk as a dense matrix 161786d161a7SShri Abhyankar */ 161886d161a7SShri Abhyankar if (nz == MATRIX_BINARY_FORMAT_DENSE) { 16199d36ed5fSBarry Smith ierr = MatLoad_MPIDense_DenseInFile(comm,fd,M,N,newmat);CHKERRQ(ierr); 162086d161a7SShri Abhyankar PetscFunctionReturn(0); 162186d161a7SShri Abhyankar } 162286d161a7SShri Abhyankar 162386d161a7SShri Abhyankar /* determine ownership of all rows */ 16242205254eSKarl Rupp if (newmat->rmap->n < 0) { 16252205254eSKarl Rupp ierr = PetscMPIIntCast(M/size + ((M % size) > rank),&m);CHKERRQ(ierr); 16262205254eSKarl Rupp } else { 16272205254eSKarl Rupp ierr = PetscMPIIntCast(newmat->rmap->n,&m);CHKERRQ(ierr); 16282205254eSKarl Rupp } 162949467e7dSBarry Smith if (newmat->cmap->n < 0) { 163049467e7dSBarry Smith n = PETSC_DECIDE; 163149467e7dSBarry Smith } else { 163249467e7dSBarry Smith ierr = PetscMPIIntCast(newmat->cmap->n,&n);CHKERRQ(ierr); 163349467e7dSBarry Smith } 163449467e7dSBarry Smith 1635854ce69bSBarry Smith ierr = PetscMalloc1(size+2,&rowners);CHKERRQ(ierr); 163686d161a7SShri Abhyankar ierr = MPI_Allgather(&m,1,MPI_INT,rowners+1,1,MPI_INT,comm);CHKERRQ(ierr); 163786d161a7SShri Abhyankar rowners[0] = 0; 163886d161a7SShri Abhyankar for (i=2; i<=size; i++) { 163986d161a7SShri Abhyankar rowners[i] += rowners[i-1]; 164086d161a7SShri Abhyankar } 164186d161a7SShri Abhyankar rstart = rowners[rank]; 164286d161a7SShri Abhyankar rend = rowners[rank+1]; 164386d161a7SShri Abhyankar 164486d161a7SShri Abhyankar /* distribute row lengths to all processors */ 164549467e7dSBarry Smith ierr = PetscMalloc1(rend-rstart,&ourlens);CHKERRQ(ierr); 164686d161a7SShri Abhyankar if (!rank) { 1647785e854fSJed Brown ierr = PetscMalloc1(M,&rowlengths);CHKERRQ(ierr); 164886d161a7SShri Abhyankar ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr); 1649785e854fSJed Brown ierr = PetscMalloc1(size,&sndcounts);CHKERRQ(ierr); 165086d161a7SShri Abhyankar for (i=0; i<size; i++) sndcounts[i] = rowners[i+1] - rowners[i]; 165186d161a7SShri Abhyankar ierr = MPI_Scatterv(rowlengths,sndcounts,rowners,MPIU_INT,ourlens,rend-rstart,MPIU_INT,0,comm);CHKERRQ(ierr); 165286d161a7SShri Abhyankar ierr = PetscFree(sndcounts);CHKERRQ(ierr); 165386d161a7SShri Abhyankar } else { 165486d161a7SShri Abhyankar ierr = MPI_Scatterv(0,0,0,MPIU_INT,ourlens,rend-rstart,MPIU_INT,0,comm);CHKERRQ(ierr); 165586d161a7SShri Abhyankar } 165686d161a7SShri Abhyankar 165786d161a7SShri Abhyankar if (!rank) { 165886d161a7SShri Abhyankar /* calculate the number of nonzeros on each processor */ 1659785e854fSJed Brown ierr = PetscMalloc1(size,&procsnz);CHKERRQ(ierr); 166086d161a7SShri Abhyankar ierr = PetscMemzero(procsnz,size*sizeof(PetscInt));CHKERRQ(ierr); 166186d161a7SShri Abhyankar for (i=0; i<size; i++) { 166286d161a7SShri Abhyankar for (j=rowners[i]; j< rowners[i+1]; j++) { 166386d161a7SShri Abhyankar procsnz[i] += rowlengths[j]; 166486d161a7SShri Abhyankar } 166586d161a7SShri Abhyankar } 166686d161a7SShri Abhyankar ierr = PetscFree(rowlengths);CHKERRQ(ierr); 166786d161a7SShri Abhyankar 166886d161a7SShri Abhyankar /* determine max buffer needed and allocate it */ 166986d161a7SShri Abhyankar maxnz = 0; 167086d161a7SShri Abhyankar for (i=0; i<size; i++) { 167186d161a7SShri Abhyankar maxnz = PetscMax(maxnz,procsnz[i]); 167286d161a7SShri Abhyankar } 1673785e854fSJed Brown ierr = PetscMalloc1(maxnz,&cols);CHKERRQ(ierr); 167486d161a7SShri Abhyankar 167586d161a7SShri Abhyankar /* read in my part of the matrix column indices */ 167686d161a7SShri Abhyankar nz = procsnz[0]; 1677785e854fSJed Brown ierr = PetscMalloc1(nz,&mycols);CHKERRQ(ierr); 167886d161a7SShri Abhyankar ierr = PetscBinaryRead(fd,mycols,nz,PETSC_INT);CHKERRQ(ierr); 167986d161a7SShri Abhyankar 168086d161a7SShri Abhyankar /* read in every one elses and ship off */ 168186d161a7SShri Abhyankar for (i=1; i<size; i++) { 168286d161a7SShri Abhyankar nz = procsnz[i]; 168386d161a7SShri Abhyankar ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT);CHKERRQ(ierr); 168486d161a7SShri Abhyankar ierr = MPI_Send(cols,nz,MPIU_INT,i,tag,comm);CHKERRQ(ierr); 168586d161a7SShri Abhyankar } 168686d161a7SShri Abhyankar ierr = PetscFree(cols);CHKERRQ(ierr); 168786d161a7SShri Abhyankar } else { 168886d161a7SShri Abhyankar /* determine buffer space needed for message */ 168986d161a7SShri Abhyankar nz = 0; 169086d161a7SShri Abhyankar for (i=0; i<m; i++) { 169186d161a7SShri Abhyankar nz += ourlens[i]; 169286d161a7SShri Abhyankar } 1693854ce69bSBarry Smith ierr = PetscMalloc1(nz+1,&mycols);CHKERRQ(ierr); 169486d161a7SShri Abhyankar 169586d161a7SShri Abhyankar /* receive message of column indices*/ 169686d161a7SShri Abhyankar ierr = MPI_Recv(mycols,nz,MPIU_INT,0,tag,comm,&status);CHKERRQ(ierr); 169786d161a7SShri Abhyankar ierr = MPI_Get_count(&status,MPIU_INT,&maxnz);CHKERRQ(ierr); 169886d161a7SShri Abhyankar if (maxnz != nz) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"something is wrong with file"); 169986d161a7SShri Abhyankar } 170086d161a7SShri Abhyankar 170149467e7dSBarry Smith ierr = MatSetSizes(newmat,m,n,M,N);CHKERRQ(ierr); 1702dfdea288SBarry Smith a = (Mat_MPIDense*)newmat->data; 17032e3566b0SBarry Smith if (!a->A || !((Mat_SeqDense*)(a->A->data))->user_alloc) { 17040298fd71SBarry Smith ierr = MatMPIDenseSetPreallocation(newmat,NULL);CHKERRQ(ierr); 1705dfdea288SBarry Smith } 170686d161a7SShri Abhyankar 170786d161a7SShri Abhyankar if (!rank) { 1708785e854fSJed Brown ierr = PetscMalloc1(maxnz,&vals);CHKERRQ(ierr); 170986d161a7SShri Abhyankar 171086d161a7SShri Abhyankar /* read in my part of the matrix numerical values */ 171186d161a7SShri Abhyankar nz = procsnz[0]; 171286d161a7SShri Abhyankar ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr); 171386d161a7SShri Abhyankar 171486d161a7SShri Abhyankar /* insert into matrix */ 171586d161a7SShri Abhyankar jj = rstart; 171686d161a7SShri Abhyankar smycols = mycols; 171786d161a7SShri Abhyankar svals = vals; 171886d161a7SShri Abhyankar for (i=0; i<m; i++) { 171986d161a7SShri Abhyankar ierr = MatSetValues(newmat,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr); 172086d161a7SShri Abhyankar smycols += ourlens[i]; 172186d161a7SShri Abhyankar svals += ourlens[i]; 172286d161a7SShri Abhyankar jj++; 172386d161a7SShri Abhyankar } 172486d161a7SShri Abhyankar 172586d161a7SShri Abhyankar /* read in other processors and ship out */ 172686d161a7SShri Abhyankar for (i=1; i<size; i++) { 172786d161a7SShri Abhyankar nz = procsnz[i]; 172886d161a7SShri Abhyankar ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr); 172986d161a7SShri Abhyankar ierr = MPI_Send(vals,nz,MPIU_SCALAR,i,((PetscObject)newmat)->tag,comm);CHKERRQ(ierr); 173086d161a7SShri Abhyankar } 173186d161a7SShri Abhyankar ierr = PetscFree(procsnz);CHKERRQ(ierr); 173286d161a7SShri Abhyankar } else { 173386d161a7SShri Abhyankar /* receive numeric values */ 1734854ce69bSBarry Smith ierr = PetscMalloc1(nz+1,&vals);CHKERRQ(ierr); 173586d161a7SShri Abhyankar 173686d161a7SShri Abhyankar /* receive message of values*/ 173786d161a7SShri Abhyankar ierr = MPI_Recv(vals,nz,MPIU_SCALAR,0,((PetscObject)newmat)->tag,comm,&status);CHKERRQ(ierr); 173886d161a7SShri Abhyankar ierr = MPI_Get_count(&status,MPIU_SCALAR,&maxnz);CHKERRQ(ierr); 173986d161a7SShri Abhyankar if (maxnz != nz) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"something is wrong with file"); 174086d161a7SShri Abhyankar 174186d161a7SShri Abhyankar /* insert into matrix */ 174286d161a7SShri Abhyankar jj = rstart; 174386d161a7SShri Abhyankar smycols = mycols; 174486d161a7SShri Abhyankar svals = vals; 174586d161a7SShri Abhyankar for (i=0; i<m; i++) { 174686d161a7SShri Abhyankar ierr = MatSetValues(newmat,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr); 174786d161a7SShri Abhyankar smycols += ourlens[i]; 174886d161a7SShri Abhyankar svals += ourlens[i]; 174986d161a7SShri Abhyankar jj++; 175086d161a7SShri Abhyankar } 175186d161a7SShri Abhyankar } 175249467e7dSBarry Smith ierr = PetscFree(ourlens);CHKERRQ(ierr); 175386d161a7SShri Abhyankar ierr = PetscFree(vals);CHKERRQ(ierr); 175486d161a7SShri Abhyankar ierr = PetscFree(mycols);CHKERRQ(ierr); 175586d161a7SShri Abhyankar ierr = PetscFree(rowners);CHKERRQ(ierr); 175686d161a7SShri Abhyankar 175786d161a7SShri Abhyankar ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 175886d161a7SShri Abhyankar ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 175986d161a7SShri Abhyankar PetscFunctionReturn(0); 176086d161a7SShri Abhyankar } 176186d161a7SShri Abhyankar 1762ace3abfcSBarry Smith PetscErrorCode MatEqual_MPIDense(Mat A,Mat B,PetscBool *flag) 17636e4ee0c6SHong Zhang { 17646e4ee0c6SHong Zhang Mat_MPIDense *matB = (Mat_MPIDense*)B->data,*matA = (Mat_MPIDense*)A->data; 17656e4ee0c6SHong Zhang Mat a,b; 1766ace3abfcSBarry Smith PetscBool flg; 17676e4ee0c6SHong Zhang PetscErrorCode ierr; 176890ace30eSBarry Smith 17696e4ee0c6SHong Zhang PetscFunctionBegin; 17706e4ee0c6SHong Zhang a = matA->A; 17716e4ee0c6SHong Zhang b = matB->A; 17726e4ee0c6SHong Zhang ierr = MatEqual(a,b,&flg);CHKERRQ(ierr); 1773b2566f29SBarry Smith ierr = MPIU_Allreduce(&flg,flag,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 17746e4ee0c6SHong Zhang PetscFunctionReturn(0); 17756e4ee0c6SHong Zhang } 177690ace30eSBarry Smith 1777baa3c1c6SHong Zhang PetscErrorCode MatDestroy_MatTransMatMult_MPIDense_MPIDense(Mat A) 1778baa3c1c6SHong Zhang { 1779baa3c1c6SHong Zhang PetscErrorCode ierr; 1780baa3c1c6SHong Zhang Mat_MPIDense *a = (Mat_MPIDense*)A->data; 1781baa3c1c6SHong Zhang Mat_TransMatMultDense *atb = a->atbdense; 1782baa3c1c6SHong Zhang 1783baa3c1c6SHong Zhang PetscFunctionBegin; 1784c5ef1628SHong Zhang ierr = PetscFree3(atb->sendbuf,atb->atbarray,atb->recvcounts);CHKERRQ(ierr); 1785baa3c1c6SHong Zhang ierr = (atb->destroy)(A);CHKERRQ(ierr); 1786baa3c1c6SHong Zhang ierr = PetscFree(atb);CHKERRQ(ierr); 1787baa3c1c6SHong Zhang PetscFunctionReturn(0); 1788baa3c1c6SHong Zhang } 1789baa3c1c6SHong Zhang 1790cb20be35SHong Zhang PetscErrorCode MatTransposeMatMultNumeric_MPIDense_MPIDense(Mat A,Mat B,Mat C) 1791cb20be35SHong Zhang { 1792baa3c1c6SHong Zhang Mat_MPIDense *a=(Mat_MPIDense*)A->data, *b=(Mat_MPIDense*)B->data, *c=(Mat_MPIDense*)C->data; 1793660d5466SHong Zhang Mat_SeqDense *aseq=(Mat_SeqDense*)(a->A)->data, *bseq=(Mat_SeqDense*)(b->A)->data; 1794baa3c1c6SHong Zhang Mat_TransMatMultDense *atb = c->atbdense; 1795cb20be35SHong Zhang PetscErrorCode ierr; 1796cb20be35SHong Zhang MPI_Comm comm; 1797d5017740SHong Zhang PetscMPIInt rank,size,*recvcounts=atb->recvcounts; 1798c5ef1628SHong Zhang PetscScalar *carray,*atbarray=atb->atbarray,*sendbuf=atb->sendbuf; 1799d5017740SHong Zhang PetscInt i,cN=C->cmap->N,cM=C->rmap->N,proc,k,j; 1800660d5466SHong Zhang PetscScalar _DOne=1.0,_DZero=0.0; 1801adc7a786SHong Zhang PetscBLASInt am,an,bn,aN; 1802e68c0b26SHong Zhang const PetscInt *ranges; 1803cb20be35SHong Zhang 1804cb20be35SHong Zhang PetscFunctionBegin; 1805cb20be35SHong Zhang ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 1806cb20be35SHong Zhang ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 1807cb20be35SHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1808e68c0b26SHong Zhang 1809c5ef1628SHong Zhang /* compute atbarray = aseq^T * bseq */ 1810660d5466SHong Zhang ierr = PetscBLASIntCast(a->A->cmap->n,&an);CHKERRQ(ierr); 1811660d5466SHong Zhang ierr = PetscBLASIntCast(b->A->cmap->n,&bn);CHKERRQ(ierr); 1812660d5466SHong Zhang ierr = PetscBLASIntCast(a->A->rmap->n,&am);CHKERRQ(ierr); 1813adc7a786SHong Zhang ierr = PetscBLASIntCast(A->cmap->N,&aN);CHKERRQ(ierr); 1814adc7a786SHong Zhang PetscStackCallBLAS("BLASgemm",BLASgemm_("T","N",&an,&bn,&am,&_DOne,aseq->v,&aseq->lda,bseq->v,&bseq->lda,&_DZero,atbarray,&aN)); 1815cb20be35SHong Zhang 1816cb20be35SHong Zhang ierr = MatGetOwnershipRanges(C,&ranges);CHKERRQ(ierr); 1817c5ef1628SHong Zhang for (i=0; i<size; i++) recvcounts[i] = (ranges[i+1] - ranges[i])*cN; 1818cb20be35SHong Zhang 1819660d5466SHong Zhang /* arrange atbarray into sendbuf */ 1820cb20be35SHong Zhang k = 0; 1821cb20be35SHong Zhang for (proc=0; proc<size; proc++) { 1822baa3c1c6SHong Zhang for (j=0; j<cN; j++) { 1823c5ef1628SHong Zhang for (i=ranges[proc]; i<ranges[proc+1]; i++) sendbuf[k++] = atbarray[i+j*cM]; 1824cb20be35SHong Zhang } 1825cb20be35SHong Zhang } 1826c5ef1628SHong Zhang /* sum all atbarray to local values of C */ 1827660d5466SHong Zhang ierr = MatDenseGetArray(c->A,&carray);CHKERRQ(ierr); 18283462b7efSHong Zhang ierr = MPI_Reduce_scatter(sendbuf,carray,recvcounts,MPIU_SCALAR,MPIU_SUM,comm);CHKERRQ(ierr); 1829660d5466SHong Zhang ierr = MatDenseRestoreArray(c->A,&carray);CHKERRQ(ierr); 1830cb20be35SHong Zhang PetscFunctionReturn(0); 1831cb20be35SHong Zhang } 1832cb20be35SHong Zhang 1833cb20be35SHong Zhang PetscErrorCode MatTransposeMatMultSymbolic_MPIDense_MPIDense(Mat A,Mat B,PetscReal fill,Mat *C) 1834cb20be35SHong Zhang { 1835cb20be35SHong Zhang PetscErrorCode ierr; 1836baa3c1c6SHong Zhang Mat Cdense; 1837cb20be35SHong Zhang MPI_Comm comm; 1838baa3c1c6SHong Zhang PetscMPIInt size; 1839660d5466SHong Zhang PetscInt cm=A->cmap->n,cM,cN=B->cmap->N; 1840baa3c1c6SHong Zhang Mat_MPIDense *c; 1841baa3c1c6SHong Zhang Mat_TransMatMultDense *atb; 1842cb20be35SHong Zhang 1843cb20be35SHong Zhang PetscFunctionBegin; 1844baa3c1c6SHong Zhang ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 1845cb20be35SHong Zhang if (A->rmap->rstart != B->rmap->rstart || A->rmap->rend != B->rmap->rend) { 1846cb20be35SHong 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); 1847cb20be35SHong Zhang } 1848cb20be35SHong Zhang 1849baa3c1c6SHong Zhang /* create matrix product Cdense */ 1850baa3c1c6SHong Zhang ierr = MatCreate(comm,&Cdense);CHKERRQ(ierr); 1851660d5466SHong Zhang ierr = MatSetSizes(Cdense,cm,B->cmap->n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr); 1852baa3c1c6SHong Zhang ierr = MatSetType(Cdense,MATMPIDENSE);CHKERRQ(ierr); 1853baa3c1c6SHong Zhang ierr = MatMPIDenseSetPreallocation(Cdense,NULL);CHKERRQ(ierr); 1854baa3c1c6SHong Zhang ierr = MatAssemblyBegin(Cdense,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1855baa3c1c6SHong Zhang ierr = MatAssemblyEnd(Cdense,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1856baa3c1c6SHong Zhang *C = Cdense; 1857baa3c1c6SHong Zhang 1858baa3c1c6SHong Zhang /* create data structure for reuse Cdense */ 1859baa3c1c6SHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1860baa3c1c6SHong Zhang ierr = PetscNew(&atb);CHKERRQ(ierr); 1861660d5466SHong Zhang cM = Cdense->rmap->N; 1862c5ef1628SHong Zhang ierr = PetscMalloc3(cM*cN,&atb->sendbuf,cM*cN,&atb->atbarray,size,&atb->recvcounts);CHKERRQ(ierr); 1863baa3c1c6SHong Zhang 1864baa3c1c6SHong Zhang c = (Mat_MPIDense*)Cdense->data; 1865baa3c1c6SHong Zhang c->atbdense = atb; 1866baa3c1c6SHong Zhang atb->destroy = Cdense->ops->destroy; 1867baa3c1c6SHong Zhang Cdense->ops->destroy = MatDestroy_MatTransMatMult_MPIDense_MPIDense; 1868cb20be35SHong Zhang PetscFunctionReturn(0); 1869cb20be35SHong Zhang } 1870cb20be35SHong Zhang 1871cb20be35SHong Zhang PetscErrorCode MatTransposeMatMult_MPIDense_MPIDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 1872cb20be35SHong Zhang { 1873cb20be35SHong Zhang PetscErrorCode ierr; 1874cb20be35SHong Zhang 1875cb20be35SHong Zhang PetscFunctionBegin; 1876cb20be35SHong Zhang if (scall == MAT_INITIAL_MATRIX) { 1877cb20be35SHong Zhang ierr = MatTransposeMatMultSymbolic_MPIDense_MPIDense(A,B,fill,C);CHKERRQ(ierr); 1878cb20be35SHong Zhang } 1879cb20be35SHong Zhang ierr = MatTransposeMatMultNumeric_MPIDense_MPIDense(A,B,*C);CHKERRQ(ierr); 1880cb20be35SHong Zhang PetscFunctionReturn(0); 1881cb20be35SHong Zhang } 1882320f2790SHong Zhang 1883320f2790SHong Zhang PetscErrorCode MatDestroy_MatMatMult_MPIDense_MPIDense(Mat A) 1884320f2790SHong Zhang { 1885320f2790SHong Zhang PetscErrorCode ierr; 1886320f2790SHong Zhang Mat_MPIDense *a = (Mat_MPIDense*)A->data; 1887320f2790SHong Zhang Mat_MatMultDense *ab = a->abdense; 1888320f2790SHong Zhang 1889320f2790SHong Zhang PetscFunctionBegin; 1890320f2790SHong Zhang ierr = MatDestroy(&ab->Ce);CHKERRQ(ierr); 1891320f2790SHong Zhang ierr = MatDestroy(&ab->Ae);CHKERRQ(ierr); 1892320f2790SHong Zhang ierr = MatDestroy(&ab->Be);CHKERRQ(ierr); 1893320f2790SHong Zhang 1894320f2790SHong Zhang ierr = (ab->destroy)(A);CHKERRQ(ierr); 1895320f2790SHong Zhang ierr = PetscFree(ab);CHKERRQ(ierr); 1896320f2790SHong Zhang PetscFunctionReturn(0); 1897320f2790SHong Zhang } 1898320f2790SHong Zhang 18995242a7b1SHong Zhang #if defined(PETSC_HAVE_ELEMENTAL) 1900320f2790SHong Zhang PetscErrorCode MatMatMultNumeric_MPIDense_MPIDense(Mat A,Mat B,Mat C) 1901320f2790SHong Zhang { 1902320f2790SHong Zhang PetscErrorCode ierr; 1903320f2790SHong Zhang Mat_MPIDense *c=(Mat_MPIDense*)C->data; 1904320f2790SHong Zhang Mat_MatMultDense *ab=c->abdense; 1905320f2790SHong Zhang 1906320f2790SHong Zhang PetscFunctionBegin; 1907de0a22f0SHong Zhang ierr = MatConvert_MPIDense_Elemental(A,MATELEMENTAL,MAT_REUSE_MATRIX, &ab->Ae);CHKERRQ(ierr); 1908de0a22f0SHong Zhang ierr = MatConvert_MPIDense_Elemental(B,MATELEMENTAL,MAT_REUSE_MATRIX, &ab->Be);CHKERRQ(ierr); 1909320f2790SHong Zhang ierr = MatMatMultNumeric(ab->Ae,ab->Be,ab->Ce);CHKERRQ(ierr); 1910de0a22f0SHong Zhang ierr = MatConvert(ab->Ce,MATMPIDENSE,MAT_REUSE_MATRIX,&C);CHKERRQ(ierr); 1911320f2790SHong Zhang PetscFunctionReturn(0); 1912320f2790SHong Zhang } 1913320f2790SHong Zhang 1914320f2790SHong Zhang PetscErrorCode MatMatMultSymbolic_MPIDense_MPIDense(Mat A,Mat B,PetscReal fill,Mat *C) 1915320f2790SHong Zhang { 1916320f2790SHong Zhang PetscErrorCode ierr; 1917320f2790SHong Zhang Mat Ae,Be,Ce; 1918320f2790SHong Zhang Mat_MPIDense *c; 1919320f2790SHong Zhang Mat_MatMultDense *ab; 1920320f2790SHong Zhang 1921320f2790SHong Zhang PetscFunctionBegin; 1922320f2790SHong Zhang if (A->cmap->rstart != B->rmap->rstart || A->cmap->rend != B->rmap->rend) { 1923378336b6SHong 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); 1924320f2790SHong Zhang } 1925320f2790SHong Zhang 1926320f2790SHong Zhang /* convert A and B to Elemental matrices Ae and Be */ 1927320f2790SHong Zhang ierr = MatConvert(A,MATELEMENTAL,MAT_INITIAL_MATRIX, &Ae);CHKERRQ(ierr); 1928320f2790SHong Zhang ierr = MatConvert(B,MATELEMENTAL,MAT_INITIAL_MATRIX, &Be);CHKERRQ(ierr); 1929320f2790SHong Zhang 1930320f2790SHong Zhang /* Ce = Ae*Be */ 1931320f2790SHong Zhang ierr = MatMatMultSymbolic(Ae,Be,fill,&Ce);CHKERRQ(ierr); 1932320f2790SHong Zhang ierr = MatMatMultNumeric(Ae,Be,Ce);CHKERRQ(ierr); 1933320f2790SHong Zhang 1934320f2790SHong Zhang /* convert Ce to C */ 1935320f2790SHong Zhang ierr = MatConvert(Ce,MATMPIDENSE,MAT_INITIAL_MATRIX,C);CHKERRQ(ierr); 1936320f2790SHong Zhang 1937320f2790SHong Zhang /* create data structure for reuse Cdense */ 1938320f2790SHong Zhang ierr = PetscNew(&ab);CHKERRQ(ierr); 1939320f2790SHong Zhang c = (Mat_MPIDense*)(*C)->data; 1940320f2790SHong Zhang c->abdense = ab; 1941320f2790SHong Zhang 1942320f2790SHong Zhang ab->Ae = Ae; 1943320f2790SHong Zhang ab->Be = Be; 1944320f2790SHong Zhang ab->Ce = Ce; 1945320f2790SHong Zhang ab->destroy = (*C)->ops->destroy; 1946320f2790SHong Zhang (*C)->ops->destroy = MatDestroy_MatMatMult_MPIDense_MPIDense; 19479775878aSHong Zhang (*C)->ops->matmultnumeric = MatMatMultNumeric_MPIDense_MPIDense; 1948320f2790SHong Zhang PetscFunctionReturn(0); 1949320f2790SHong Zhang } 1950320f2790SHong Zhang 1951150d2497SBarry Smith PETSC_INTERN PetscErrorCode MatMatMult_MPIDense_MPIDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 1952320f2790SHong Zhang { 1953320f2790SHong Zhang PetscErrorCode ierr; 1954320f2790SHong Zhang 1955320f2790SHong Zhang PetscFunctionBegin; 1956320f2790SHong Zhang if (scall == MAT_INITIAL_MATRIX) { /* simbolic product includes numeric product */ 195757299f2fSHong Zhang ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr); 1958320f2790SHong Zhang ierr = MatMatMultSymbolic_MPIDense_MPIDense(A,B,fill,C);CHKERRQ(ierr); 195957299f2fSHong Zhang ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr); 1960320f2790SHong Zhang } else { 196157299f2fSHong Zhang ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr); 1962320f2790SHong Zhang ierr = MatMatMultNumeric_MPIDense_MPIDense(A,B,*C);CHKERRQ(ierr); 196357299f2fSHong Zhang ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr); 1964320f2790SHong Zhang } 1965320f2790SHong Zhang PetscFunctionReturn(0); 1966320f2790SHong Zhang } 19675242a7b1SHong Zhang #endif 1968