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 14749a6ff4bSBarry Smith static PetscErrorCode MatDenseGetLDA_MPIDense(Mat A,PetscInt *lda) 14849a6ff4bSBarry Smith { 14949a6ff4bSBarry Smith Mat_MPIDense *a = (Mat_MPIDense*)A->data; 15049a6ff4bSBarry Smith PetscErrorCode ierr; 15149a6ff4bSBarry Smith 15249a6ff4bSBarry Smith PetscFunctionBegin; 15349a6ff4bSBarry Smith ierr = MatDenseGetLDA(a->A,lda);CHKERRQ(ierr); 15449a6ff4bSBarry Smith PetscFunctionReturn(0); 15549a6ff4bSBarry Smith } 15649a6ff4bSBarry Smith 157d3042a70SBarry Smith static PetscErrorCode MatDenseGetArray_MPIDense(Mat A,PetscScalar *array[]) 158ff14e315SSatish Balay { 159ff14e315SSatish Balay Mat_MPIDense *a = (Mat_MPIDense*)A->data; 160dfbe8321SBarry Smith PetscErrorCode ierr; 161ff14e315SSatish Balay 1623a40ed3dSBarry Smith PetscFunctionBegin; 1638c778c55SBarry Smith ierr = MatDenseGetArray(a->A,array);CHKERRQ(ierr); 1643a40ed3dSBarry Smith PetscFunctionReturn(0); 165ff14e315SSatish Balay } 166ff14e315SSatish Balay 1678572280aSBarry Smith static PetscErrorCode MatDenseGetArrayRead_MPIDense(Mat A,const PetscScalar *array[]) 1688572280aSBarry Smith { 1698572280aSBarry Smith Mat_MPIDense *a = (Mat_MPIDense*)A->data; 1708572280aSBarry Smith PetscErrorCode ierr; 1718572280aSBarry Smith 1728572280aSBarry Smith PetscFunctionBegin; 1738572280aSBarry Smith ierr = MatDenseGetArrayRead(a->A,array);CHKERRQ(ierr); 1748572280aSBarry Smith PetscFunctionReturn(0); 1758572280aSBarry Smith } 1768572280aSBarry Smith 177d3042a70SBarry Smith static PetscErrorCode MatDensePlaceArray_MPIDense(Mat A,const PetscScalar array[]) 178d3042a70SBarry Smith { 179d3042a70SBarry Smith Mat_MPIDense *a = (Mat_MPIDense*)A->data; 180d3042a70SBarry Smith PetscErrorCode ierr; 181d3042a70SBarry Smith 182d3042a70SBarry Smith PetscFunctionBegin; 183d3042a70SBarry Smith ierr = MatDensePlaceArray(a->A,array);CHKERRQ(ierr); 184d3042a70SBarry Smith PetscFunctionReturn(0); 185d3042a70SBarry Smith } 186d3042a70SBarry Smith 187d3042a70SBarry Smith static PetscErrorCode MatDenseResetArray_MPIDense(Mat A) 188d3042a70SBarry Smith { 189d3042a70SBarry Smith Mat_MPIDense *a = (Mat_MPIDense*)A->data; 190d3042a70SBarry Smith PetscErrorCode ierr; 191d3042a70SBarry Smith 192d3042a70SBarry Smith PetscFunctionBegin; 193d3042a70SBarry Smith ierr = MatDenseResetArray(a->A);CHKERRQ(ierr); 194d3042a70SBarry Smith PetscFunctionReturn(0); 195d3042a70SBarry Smith } 196d3042a70SBarry Smith 1977dae84e0SHong Zhang static PetscErrorCode MatCreateSubMatrix_MPIDense(Mat A,IS isrow,IS iscol,MatReuse scall,Mat *B) 198ca3fa75bSLois Curfman McInnes { 199ca3fa75bSLois Curfman McInnes Mat_MPIDense *mat = (Mat_MPIDense*)A->data,*newmatd; 200ca3fa75bSLois Curfman McInnes Mat_SeqDense *lmat = (Mat_SeqDense*)mat->A->data; 2016849ba73SBarry Smith PetscErrorCode ierr; 2024aa3045dSJed Brown PetscInt i,j,rstart,rend,nrows,ncols,Ncols,nlrows,nlcols; 2035d0c19d7SBarry Smith const PetscInt *irow,*icol; 20487828ca2SBarry Smith PetscScalar *av,*bv,*v = lmat->v; 205ca3fa75bSLois Curfman McInnes Mat newmat; 2064aa3045dSJed Brown IS iscol_local; 20742a884f0SBarry Smith MPI_Comm comm_is,comm_mat; 208ca3fa75bSLois Curfman McInnes 209ca3fa75bSLois Curfman McInnes PetscFunctionBegin; 21042a884f0SBarry Smith ierr = PetscObjectGetComm((PetscObject)A,&comm_mat);CHKERRQ(ierr); 21142a884f0SBarry Smith ierr = PetscObjectGetComm((PetscObject)iscol,&comm_is);CHKERRQ(ierr); 21242a884f0SBarry Smith if (comm_mat != comm_is) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NOTSAMECOMM,"IS communicator must match matrix communicator"); 21342a884f0SBarry Smith 2144aa3045dSJed Brown ierr = ISAllGather(iscol,&iscol_local);CHKERRQ(ierr); 215ca3fa75bSLois Curfman McInnes ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr); 2164aa3045dSJed Brown ierr = ISGetIndices(iscol_local,&icol);CHKERRQ(ierr); 217b9b97703SBarry Smith ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr); 218b9b97703SBarry Smith ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr); 2194aa3045dSJed Brown ierr = ISGetSize(iscol,&Ncols);CHKERRQ(ierr); /* global number of columns, size of iscol_local */ 220ca3fa75bSLois Curfman McInnes 221ca3fa75bSLois Curfman McInnes /* No parallel redistribution currently supported! Should really check each index set 2227eba5e9cSLois Curfman McInnes to comfirm that it is OK. ... Currently supports only submatrix same partitioning as 2237eba5e9cSLois Curfman McInnes original matrix! */ 224ca3fa75bSLois Curfman McInnes 225ca3fa75bSLois Curfman McInnes ierr = MatGetLocalSize(A,&nlrows,&nlcols);CHKERRQ(ierr); 2267eba5e9cSLois Curfman McInnes ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr); 227ca3fa75bSLois Curfman McInnes 228ca3fa75bSLois Curfman McInnes /* Check submatrix call */ 229ca3fa75bSLois Curfman McInnes if (scall == MAT_REUSE_MATRIX) { 230e32f2f54SBarry Smith /* SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Reused submatrix wrong size"); */ 2317eba5e9cSLois Curfman McInnes /* Really need to test rows and column sizes! */ 232ca3fa75bSLois Curfman McInnes newmat = *B; 233ca3fa75bSLois Curfman McInnes } else { 234ca3fa75bSLois Curfman McInnes /* Create and fill new matrix */ 235ce94432eSBarry Smith ierr = MatCreate(PetscObjectComm((PetscObject)A),&newmat);CHKERRQ(ierr); 2364aa3045dSJed Brown ierr = MatSetSizes(newmat,nrows,ncols,PETSC_DECIDE,Ncols);CHKERRQ(ierr); 2377adad957SLisandro Dalcin ierr = MatSetType(newmat,((PetscObject)A)->type_name);CHKERRQ(ierr); 2380298fd71SBarry Smith ierr = MatMPIDenseSetPreallocation(newmat,NULL);CHKERRQ(ierr); 239ca3fa75bSLois Curfman McInnes } 240ca3fa75bSLois Curfman McInnes 241ca3fa75bSLois Curfman McInnes /* Now extract the data pointers and do the copy, column at a time */ 242ca3fa75bSLois Curfman McInnes newmatd = (Mat_MPIDense*)newmat->data; 243ca3fa75bSLois Curfman McInnes bv = ((Mat_SeqDense*)newmatd->A->data)->v; 244ca3fa75bSLois Curfman McInnes 2454aa3045dSJed Brown for (i=0; i<Ncols; i++) { 24625a33276SHong Zhang av = v + ((Mat_SeqDense*)mat->A->data)->lda*icol[i]; 247ca3fa75bSLois Curfman McInnes for (j=0; j<nrows; j++) { 2487eba5e9cSLois Curfman McInnes *bv++ = av[irow[j] - rstart]; 249ca3fa75bSLois Curfman McInnes } 250ca3fa75bSLois Curfman McInnes } 251ca3fa75bSLois Curfman McInnes 252ca3fa75bSLois Curfman McInnes /* Assemble the matrices so that the correct flags are set */ 253ca3fa75bSLois Curfman McInnes ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 254ca3fa75bSLois Curfman McInnes ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 255ca3fa75bSLois Curfman McInnes 256ca3fa75bSLois Curfman McInnes /* Free work space */ 257ca3fa75bSLois Curfman McInnes ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr); 2585bdf786aSShri Abhyankar ierr = ISRestoreIndices(iscol_local,&icol);CHKERRQ(ierr); 25932bb1f2dSLisandro Dalcin ierr = ISDestroy(&iscol_local);CHKERRQ(ierr); 260ca3fa75bSLois Curfman McInnes *B = newmat; 261ca3fa75bSLois Curfman McInnes PetscFunctionReturn(0); 262ca3fa75bSLois Curfman McInnes } 263ca3fa75bSLois Curfman McInnes 2648c778c55SBarry Smith PetscErrorCode MatDenseRestoreArray_MPIDense(Mat A,PetscScalar *array[]) 265ff14e315SSatish Balay { 26673a71a0fSBarry Smith Mat_MPIDense *a = (Mat_MPIDense*)A->data; 26773a71a0fSBarry Smith PetscErrorCode ierr; 26873a71a0fSBarry Smith 2693a40ed3dSBarry Smith PetscFunctionBegin; 2708c778c55SBarry Smith ierr = MatDenseRestoreArray(a->A,array);CHKERRQ(ierr); 2713a40ed3dSBarry Smith PetscFunctionReturn(0); 272ff14e315SSatish Balay } 273ff14e315SSatish Balay 2748572280aSBarry Smith PetscErrorCode MatDenseRestoreArrayRead_MPIDense(Mat A,const PetscScalar *array[]) 2758572280aSBarry Smith { 2768572280aSBarry Smith Mat_MPIDense *a = (Mat_MPIDense*)A->data; 2778572280aSBarry Smith PetscErrorCode ierr; 2788572280aSBarry Smith 2798572280aSBarry Smith PetscFunctionBegin; 2808572280aSBarry Smith ierr = MatDenseRestoreArrayRead(a->A,array);CHKERRQ(ierr); 2818572280aSBarry Smith PetscFunctionReturn(0); 2828572280aSBarry Smith } 2838572280aSBarry Smith 284dfbe8321SBarry Smith PetscErrorCode MatAssemblyBegin_MPIDense(Mat mat,MatAssemblyType mode) 2858965ea79SLois Curfman McInnes { 28639ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 287ce94432eSBarry Smith MPI_Comm comm; 288dfbe8321SBarry Smith PetscErrorCode ierr; 28913f74950SBarry Smith PetscInt nstash,reallocs; 2908965ea79SLois Curfman McInnes InsertMode addv; 2918965ea79SLois Curfman McInnes 2923a40ed3dSBarry Smith PetscFunctionBegin; 293ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)mat,&comm);CHKERRQ(ierr); 2948965ea79SLois Curfman McInnes /* make sure all processors are either in INSERTMODE or ADDMODE */ 295b2566f29SBarry Smith ierr = MPIU_Allreduce((PetscEnum*)&mat->insertmode,(PetscEnum*)&addv,1,MPIU_ENUM,MPI_BOR,comm);CHKERRQ(ierr); 296ce94432eSBarry Smith if (addv == (ADD_VALUES|INSERT_VALUES)) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Cannot mix adds/inserts on different procs"); 297e0fa3b82SLois Curfman McInnes mat->insertmode = addv; /* in case this processor had no cache */ 2988965ea79SLois Curfman McInnes 299d0f46423SBarry Smith ierr = MatStashScatterBegin_Private(mat,&mat->stash,mat->rmap->range);CHKERRQ(ierr); 3008798bf22SSatish Balay ierr = MatStashGetInfo_Private(&mat->stash,&nstash,&reallocs);CHKERRQ(ierr); 301ae15b995SBarry Smith ierr = PetscInfo2(mdn->A,"Stash has %D entries, uses %D mallocs.\n",nstash,reallocs);CHKERRQ(ierr); 3023a40ed3dSBarry Smith PetscFunctionReturn(0); 3038965ea79SLois Curfman McInnes } 3048965ea79SLois Curfman McInnes 305dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_MPIDense(Mat mat,MatAssemblyType mode) 3068965ea79SLois Curfman McInnes { 30739ddd567SLois Curfman McInnes Mat_MPIDense *mdn=(Mat_MPIDense*)mat->data; 3086849ba73SBarry Smith PetscErrorCode ierr; 30913f74950SBarry Smith PetscInt i,*row,*col,flg,j,rstart,ncols; 31013f74950SBarry Smith PetscMPIInt n; 31187828ca2SBarry Smith PetscScalar *val; 3128965ea79SLois Curfman McInnes 3133a40ed3dSBarry Smith PetscFunctionBegin; 3148965ea79SLois Curfman McInnes /* wait on receives */ 3157ef1d9bdSSatish Balay while (1) { 3168798bf22SSatish Balay ierr = MatStashScatterGetMesg_Private(&mat->stash,&n,&row,&col,&val,&flg);CHKERRQ(ierr); 3177ef1d9bdSSatish Balay if (!flg) break; 3188965ea79SLois Curfman McInnes 3197ef1d9bdSSatish Balay for (i=0; i<n;) { 3207ef1d9bdSSatish Balay /* Now identify the consecutive vals belonging to the same row */ 3212205254eSKarl Rupp for (j=i,rstart=row[j]; j<n; j++) { 3222205254eSKarl Rupp if (row[j] != rstart) break; 3232205254eSKarl Rupp } 3247ef1d9bdSSatish Balay if (j < n) ncols = j-i; 3257ef1d9bdSSatish Balay else ncols = n-i; 3267ef1d9bdSSatish Balay /* Now assemble all these values with a single function call */ 3274b4eb8d3SJed Brown ierr = MatSetValues_MPIDense(mat,1,row+i,ncols,col+i,val+i,mat->insertmode);CHKERRQ(ierr); 3287ef1d9bdSSatish Balay i = j; 3298965ea79SLois Curfman McInnes } 3307ef1d9bdSSatish Balay } 3318798bf22SSatish Balay ierr = MatStashScatterEnd_Private(&mat->stash);CHKERRQ(ierr); 3328965ea79SLois Curfman McInnes 33339ddd567SLois Curfman McInnes ierr = MatAssemblyBegin(mdn->A,mode);CHKERRQ(ierr); 33439ddd567SLois Curfman McInnes ierr = MatAssemblyEnd(mdn->A,mode);CHKERRQ(ierr); 3358965ea79SLois Curfman McInnes 3366d4a8577SBarry Smith if (!mat->was_assembled && mode == MAT_FINAL_ASSEMBLY) { 33739ddd567SLois Curfman McInnes ierr = MatSetUpMultiply_MPIDense(mat);CHKERRQ(ierr); 3388965ea79SLois Curfman McInnes } 3393a40ed3dSBarry Smith PetscFunctionReturn(0); 3408965ea79SLois Curfman McInnes } 3418965ea79SLois Curfman McInnes 342dfbe8321SBarry Smith PetscErrorCode MatZeroEntries_MPIDense(Mat A) 3438965ea79SLois Curfman McInnes { 344dfbe8321SBarry Smith PetscErrorCode ierr; 34539ddd567SLois Curfman McInnes Mat_MPIDense *l = (Mat_MPIDense*)A->data; 3463a40ed3dSBarry Smith 3473a40ed3dSBarry Smith PetscFunctionBegin; 3483a40ed3dSBarry Smith ierr = MatZeroEntries(l->A);CHKERRQ(ierr); 3493a40ed3dSBarry Smith PetscFunctionReturn(0); 3508965ea79SLois Curfman McInnes } 3518965ea79SLois Curfman McInnes 3528965ea79SLois Curfman McInnes /* the code does not do the diagonal entries correctly unless the 3538965ea79SLois Curfman McInnes matrix is square and the column and row owerships are identical. 3548965ea79SLois Curfman McInnes This is a BUG. The only way to fix it seems to be to access 3553501a2bdSLois Curfman McInnes mdn->A and mdn->B directly and not through the MatZeroRows() 3568965ea79SLois Curfman McInnes routine. 3578965ea79SLois Curfman McInnes */ 3582b40b63fSBarry Smith PetscErrorCode MatZeroRows_MPIDense(Mat A,PetscInt N,const PetscInt rows[],PetscScalar diag,Vec x,Vec b) 3598965ea79SLois Curfman McInnes { 36039ddd567SLois Curfman McInnes Mat_MPIDense *l = (Mat_MPIDense*)A->data; 3616849ba73SBarry Smith PetscErrorCode ierr; 362d0f46423SBarry Smith PetscInt i,*owners = A->rmap->range; 36376ec1555SBarry Smith PetscInt *sizes,j,idx,nsends; 36413f74950SBarry Smith PetscInt nmax,*svalues,*starts,*owner,nrecvs; 3657adad957SLisandro Dalcin PetscInt *rvalues,tag = ((PetscObject)A)->tag,count,base,slen,*source; 36613f74950SBarry Smith PetscInt *lens,*lrows,*values; 36713f74950SBarry Smith PetscMPIInt n,imdex,rank = l->rank,size = l->size; 368ce94432eSBarry Smith MPI_Comm comm; 3698965ea79SLois Curfman McInnes MPI_Request *send_waits,*recv_waits; 3708965ea79SLois Curfman McInnes MPI_Status recv_status,*send_status; 371ace3abfcSBarry Smith PetscBool found; 37297b48c8fSBarry Smith const PetscScalar *xx; 37397b48c8fSBarry Smith PetscScalar *bb; 3748965ea79SLois Curfman McInnes 3753a40ed3dSBarry Smith PetscFunctionBegin; 376ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 377ce94432eSBarry Smith if (A->rmap->N != A->cmap->N) SETERRQ(comm,PETSC_ERR_SUP,"Only handles square matrices"); 378b9679d65SBarry Smith if (A->rmap->n != A->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only handles matrices with identical column and row ownership"); 3798965ea79SLois Curfman McInnes /* first count number of contributors to each processor */ 380f628708eSJed Brown ierr = PetscCalloc1(2*size,&sizes);CHKERRQ(ierr); 381037dbc42SBarry Smith ierr = PetscMalloc1(N+1,&owner);CHKERRQ(ierr); /* see note*/ 3828965ea79SLois Curfman McInnes for (i=0; i<N; i++) { 3838965ea79SLois Curfman McInnes idx = rows[i]; 38435d8aa7fSBarry Smith found = PETSC_FALSE; 3858965ea79SLois Curfman McInnes for (j=0; j<size; j++) { 3868965ea79SLois Curfman McInnes if (idx >= owners[j] && idx < owners[j+1]) { 38776ec1555SBarry Smith sizes[2*j]++; sizes[2*j+1] = 1; owner[i] = j; found = PETSC_TRUE; break; 3888965ea79SLois Curfman McInnes } 3898965ea79SLois Curfman McInnes } 390e32f2f54SBarry Smith if (!found) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Index out of range"); 3918965ea79SLois Curfman McInnes } 3922205254eSKarl Rupp nsends = 0; 39376ec1555SBarry Smith for (i=0; i<size; i++) nsends += sizes[2*i+1]; 3948965ea79SLois Curfman McInnes 3958965ea79SLois Curfman McInnes /* inform other processors of number of messages and max length*/ 39676ec1555SBarry Smith ierr = PetscMaxSum(comm,sizes,&nmax,&nrecvs);CHKERRQ(ierr); 3978965ea79SLois Curfman McInnes 3988965ea79SLois Curfman McInnes /* post receives: */ 399785e854fSJed Brown ierr = PetscMalloc1((nrecvs+1)*(nmax+1),&rvalues);CHKERRQ(ierr); 400854ce69bSBarry Smith ierr = PetscMalloc1(nrecvs+1,&recv_waits);CHKERRQ(ierr); 4018965ea79SLois Curfman McInnes for (i=0; i<nrecvs; i++) { 40213f74950SBarry Smith ierr = MPI_Irecv(rvalues+nmax*i,nmax,MPIU_INT,MPI_ANY_SOURCE,tag,comm,recv_waits+i);CHKERRQ(ierr); 4038965ea79SLois Curfman McInnes } 4048965ea79SLois Curfman McInnes 4058965ea79SLois Curfman McInnes /* do sends: 4068965ea79SLois Curfman McInnes 1) starts[i] gives the starting index in svalues for stuff going to 4078965ea79SLois Curfman McInnes the ith processor 4088965ea79SLois Curfman McInnes */ 409854ce69bSBarry Smith ierr = PetscMalloc1(N+1,&svalues);CHKERRQ(ierr); 410854ce69bSBarry Smith ierr = PetscMalloc1(nsends+1,&send_waits);CHKERRQ(ierr); 411854ce69bSBarry Smith ierr = PetscMalloc1(size+1,&starts);CHKERRQ(ierr); 4128965ea79SLois Curfman McInnes 4138965ea79SLois Curfman McInnes starts[0] = 0; 41476ec1555SBarry Smith for (i=1; i<size; i++) starts[i] = starts[i-1] + sizes[2*i-2]; 4152205254eSKarl Rupp for (i=0; i<N; i++) svalues[starts[owner[i]]++] = rows[i]; 4162205254eSKarl Rupp 4172205254eSKarl Rupp starts[0] = 0; 41876ec1555SBarry Smith for (i=1; i<size+1; i++) starts[i] = starts[i-1] + sizes[2*i-2]; 4198965ea79SLois Curfman McInnes count = 0; 4208965ea79SLois Curfman McInnes for (i=0; i<size; i++) { 42176ec1555SBarry Smith if (sizes[2*i+1]) { 42276ec1555SBarry Smith ierr = MPI_Isend(svalues+starts[i],sizes[2*i],MPIU_INT,i,tag,comm,send_waits+count++);CHKERRQ(ierr); 4238965ea79SLois Curfman McInnes } 4248965ea79SLois Curfman McInnes } 425606d414cSSatish Balay ierr = PetscFree(starts);CHKERRQ(ierr); 4268965ea79SLois Curfman McInnes 4278965ea79SLois Curfman McInnes base = owners[rank]; 4288965ea79SLois Curfman McInnes 4298965ea79SLois Curfman McInnes /* wait on receives */ 430dcca6d9dSJed Brown ierr = PetscMalloc2(nrecvs,&lens,nrecvs,&source);CHKERRQ(ierr); 43174ed9c26SBarry Smith count = nrecvs; 43274ed9c26SBarry Smith slen = 0; 4338965ea79SLois Curfman McInnes while (count) { 434ca161407SBarry Smith ierr = MPI_Waitany(nrecvs,recv_waits,&imdex,&recv_status);CHKERRQ(ierr); 4358965ea79SLois Curfman McInnes /* unpack receives into our local space */ 43613f74950SBarry Smith ierr = MPI_Get_count(&recv_status,MPIU_INT,&n);CHKERRQ(ierr); 4372205254eSKarl Rupp 4388965ea79SLois Curfman McInnes source[imdex] = recv_status.MPI_SOURCE; 4398965ea79SLois Curfman McInnes lens[imdex] = n; 4408965ea79SLois Curfman McInnes slen += n; 4418965ea79SLois Curfman McInnes count--; 4428965ea79SLois Curfman McInnes } 443606d414cSSatish Balay ierr = PetscFree(recv_waits);CHKERRQ(ierr); 4448965ea79SLois Curfman McInnes 4458965ea79SLois Curfman McInnes /* move the data into the send scatter */ 446854ce69bSBarry Smith ierr = PetscMalloc1(slen+1,&lrows);CHKERRQ(ierr); 4478965ea79SLois Curfman McInnes count = 0; 4488965ea79SLois Curfman McInnes for (i=0; i<nrecvs; i++) { 4498965ea79SLois Curfman McInnes values = rvalues + i*nmax; 4508965ea79SLois Curfman McInnes for (j=0; j<lens[i]; j++) { 4518965ea79SLois Curfman McInnes lrows[count++] = values[j] - base; 4528965ea79SLois Curfman McInnes } 4538965ea79SLois Curfman McInnes } 454606d414cSSatish Balay ierr = PetscFree(rvalues);CHKERRQ(ierr); 45574ed9c26SBarry Smith ierr = PetscFree2(lens,source);CHKERRQ(ierr); 456606d414cSSatish Balay ierr = PetscFree(owner);CHKERRQ(ierr); 45776ec1555SBarry Smith ierr = PetscFree(sizes);CHKERRQ(ierr); 4588965ea79SLois Curfman McInnes 45997b48c8fSBarry Smith /* fix right hand side if needed */ 46097b48c8fSBarry Smith if (x && b) { 46197b48c8fSBarry Smith ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr); 46297b48c8fSBarry Smith ierr = VecGetArray(b,&bb);CHKERRQ(ierr); 46397b48c8fSBarry Smith for (i=0; i<slen; i++) { 46497b48c8fSBarry Smith bb[lrows[i]] = diag*xx[lrows[i]]; 46597b48c8fSBarry Smith } 46697b48c8fSBarry Smith ierr = VecRestoreArrayRead(x,&xx);CHKERRQ(ierr); 46797b48c8fSBarry Smith ierr = VecRestoreArray(b,&bb);CHKERRQ(ierr); 46897b48c8fSBarry Smith } 46997b48c8fSBarry Smith 4708965ea79SLois Curfman McInnes /* actually zap the local rows */ 471b9679d65SBarry Smith ierr = MatZeroRows(l->A,slen,lrows,0.0,0,0);CHKERRQ(ierr); 472e2eb51b1SBarry Smith if (diag != 0.0) { 473b9679d65SBarry Smith Mat_SeqDense *ll = (Mat_SeqDense*)l->A->data; 474b9679d65SBarry Smith PetscInt m = ll->lda, i; 475b9679d65SBarry Smith 476b9679d65SBarry Smith for (i=0; i<slen; i++) { 477b9679d65SBarry Smith ll->v[lrows[i] + m*(A->cmap->rstart + lrows[i])] = diag; 478b9679d65SBarry Smith } 479b9679d65SBarry Smith } 480606d414cSSatish Balay ierr = PetscFree(lrows);CHKERRQ(ierr); 4818965ea79SLois Curfman McInnes 4828965ea79SLois Curfman McInnes /* wait on sends */ 4838965ea79SLois Curfman McInnes if (nsends) { 484785e854fSJed Brown ierr = PetscMalloc1(nsends,&send_status);CHKERRQ(ierr); 485ca161407SBarry Smith ierr = MPI_Waitall(nsends,send_waits,send_status);CHKERRQ(ierr); 486606d414cSSatish Balay ierr = PetscFree(send_status);CHKERRQ(ierr); 4878965ea79SLois Curfman McInnes } 488606d414cSSatish Balay ierr = PetscFree(send_waits);CHKERRQ(ierr); 489606d414cSSatish Balay ierr = PetscFree(svalues);CHKERRQ(ierr); 4903a40ed3dSBarry Smith PetscFunctionReturn(0); 4918965ea79SLois Curfman McInnes } 4928965ea79SLois Curfman McInnes 493cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatMult_SeqDense(Mat,Vec,Vec); 494cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatMultAdd_SeqDense(Mat,Vec,Vec,Vec); 495cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatMultTranspose_SeqDense(Mat,Vec,Vec); 496cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatMultTransposeAdd_SeqDense(Mat,Vec,Vec,Vec); 497cc2e6a90SBarry Smith 498dfbe8321SBarry Smith PetscErrorCode MatMult_MPIDense(Mat mat,Vec xx,Vec yy) 4998965ea79SLois Curfman McInnes { 50039ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 501dfbe8321SBarry Smith PetscErrorCode ierr; 502c456f294SBarry Smith 5033a40ed3dSBarry Smith PetscFunctionBegin; 504ca9f406cSSatish Balay ierr = VecScatterBegin(mdn->Mvctx,xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 505ca9f406cSSatish Balay ierr = VecScatterEnd(mdn->Mvctx,xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 506cc2e6a90SBarry Smith ierr = MatMult_SeqDense(mdn->A,mdn->lvec,yy);CHKERRQ(ierr); 5073a40ed3dSBarry Smith PetscFunctionReturn(0); 5088965ea79SLois Curfman McInnes } 5098965ea79SLois Curfman McInnes 510dfbe8321SBarry Smith PetscErrorCode MatMultAdd_MPIDense(Mat mat,Vec xx,Vec yy,Vec zz) 5118965ea79SLois Curfman McInnes { 51239ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 513dfbe8321SBarry Smith PetscErrorCode ierr; 514c456f294SBarry Smith 5153a40ed3dSBarry Smith PetscFunctionBegin; 516ca9f406cSSatish Balay ierr = VecScatterBegin(mdn->Mvctx,xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 517ca9f406cSSatish Balay ierr = VecScatterEnd(mdn->Mvctx,xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 518cc2e6a90SBarry Smith ierr = MatMultAdd_SeqDense(mdn->A,mdn->lvec,yy,zz);CHKERRQ(ierr); 5193a40ed3dSBarry Smith PetscFunctionReturn(0); 5208965ea79SLois Curfman McInnes } 5218965ea79SLois Curfman McInnes 522dfbe8321SBarry Smith PetscErrorCode MatMultTranspose_MPIDense(Mat A,Vec xx,Vec yy) 523096963f5SLois Curfman McInnes { 524096963f5SLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense*)A->data; 525dfbe8321SBarry Smith PetscErrorCode ierr; 52687828ca2SBarry Smith PetscScalar zero = 0.0; 527096963f5SLois Curfman McInnes 5283a40ed3dSBarry Smith PetscFunctionBegin; 5292dcb1b2aSMatthew Knepley ierr = VecSet(yy,zero);CHKERRQ(ierr); 530cc2e6a90SBarry Smith ierr = MatMultTranspose_SeqDense(a->A,xx,a->lvec);CHKERRQ(ierr); 531ca9f406cSSatish Balay ierr = VecScatterBegin(a->Mvctx,a->lvec,yy,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 532ca9f406cSSatish Balay ierr = VecScatterEnd(a->Mvctx,a->lvec,yy,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 5333a40ed3dSBarry Smith PetscFunctionReturn(0); 534096963f5SLois Curfman McInnes } 535096963f5SLois Curfman McInnes 536dfbe8321SBarry Smith PetscErrorCode MatMultTransposeAdd_MPIDense(Mat A,Vec xx,Vec yy,Vec zz) 537096963f5SLois Curfman McInnes { 538096963f5SLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense*)A->data; 539dfbe8321SBarry Smith PetscErrorCode ierr; 540096963f5SLois Curfman McInnes 5413a40ed3dSBarry Smith PetscFunctionBegin; 5423501a2bdSLois Curfman McInnes ierr = VecCopy(yy,zz);CHKERRQ(ierr); 543cc2e6a90SBarry Smith ierr = MatMultTranspose_SeqDense(a->A,xx,a->lvec);CHKERRQ(ierr); 544ca9f406cSSatish Balay ierr = VecScatterBegin(a->Mvctx,a->lvec,zz,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 545ca9f406cSSatish Balay ierr = VecScatterEnd(a->Mvctx,a->lvec,zz,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 5463a40ed3dSBarry Smith PetscFunctionReturn(0); 547096963f5SLois Curfman McInnes } 548096963f5SLois Curfman McInnes 549dfbe8321SBarry Smith PetscErrorCode MatGetDiagonal_MPIDense(Mat A,Vec v) 5508965ea79SLois Curfman McInnes { 55139ddd567SLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense*)A->data; 552096963f5SLois Curfman McInnes Mat_SeqDense *aloc = (Mat_SeqDense*)a->A->data; 553dfbe8321SBarry Smith PetscErrorCode ierr; 554d0f46423SBarry Smith PetscInt len,i,n,m = A->rmap->n,radd; 55587828ca2SBarry Smith PetscScalar *x,zero = 0.0; 556ed3cc1f0SBarry Smith 5573a40ed3dSBarry Smith PetscFunctionBegin; 5582dcb1b2aSMatthew Knepley ierr = VecSet(v,zero);CHKERRQ(ierr); 5591ebc52fbSHong Zhang ierr = VecGetArray(v,&x);CHKERRQ(ierr); 560096963f5SLois Curfman McInnes ierr = VecGetSize(v,&n);CHKERRQ(ierr); 561e32f2f54SBarry Smith if (n != A->rmap->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec"); 562d0f46423SBarry Smith len = PetscMin(a->A->rmap->n,a->A->cmap->n); 563d0f46423SBarry Smith radd = A->rmap->rstart*m; 56444cd7ae7SLois Curfman McInnes for (i=0; i<len; i++) { 565096963f5SLois Curfman McInnes x[i] = aloc->v[radd + i*m + i]; 566096963f5SLois Curfman McInnes } 5671ebc52fbSHong Zhang ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 5683a40ed3dSBarry Smith PetscFunctionReturn(0); 5698965ea79SLois Curfman McInnes } 5708965ea79SLois Curfman McInnes 571dfbe8321SBarry Smith PetscErrorCode MatDestroy_MPIDense(Mat mat) 5728965ea79SLois Curfman McInnes { 5733501a2bdSLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 574dfbe8321SBarry Smith PetscErrorCode ierr; 575ed3cc1f0SBarry Smith 5763a40ed3dSBarry Smith PetscFunctionBegin; 577aa482453SBarry Smith #if defined(PETSC_USE_LOG) 578d0f46423SBarry Smith PetscLogObjectState((PetscObject)mat,"Rows=%D, Cols=%D",mat->rmap->N,mat->cmap->N); 5798965ea79SLois Curfman McInnes #endif 5808798bf22SSatish Balay ierr = MatStashDestroy_Private(&mat->stash);CHKERRQ(ierr); 5816bf464f9SBarry Smith ierr = MatDestroy(&mdn->A);CHKERRQ(ierr); 5826bf464f9SBarry Smith ierr = VecDestroy(&mdn->lvec);CHKERRQ(ierr); 5836bf464f9SBarry Smith ierr = VecScatterDestroy(&mdn->Mvctx);CHKERRQ(ierr); 58401b82886SBarry Smith 585bf0cc555SLisandro Dalcin ierr = PetscFree(mat->data);CHKERRQ(ierr); 586dbd8c25aSHong Zhang ierr = PetscObjectChangeTypeName((PetscObject)mat,0);CHKERRQ(ierr); 5878baccfbdSHong Zhang 58849a6ff4bSBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatDenseGetLDA_C",NULL);CHKERRQ(ierr); 5898baccfbdSHong Zhang ierr = PetscObjectComposeFunction((PetscObject)mat,"MatDenseGetArray_C",NULL);CHKERRQ(ierr); 5908572280aSBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatDenseRestoreArray_C",NULL);CHKERRQ(ierr); 5918572280aSBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatDenseGetArrayRead_C",NULL);CHKERRQ(ierr); 5928572280aSBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatDenseRestoreArrayRead_C",NULL);CHKERRQ(ierr); 593d3042a70SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatDensePlaceArray_C",NULL);CHKERRQ(ierr); 594d3042a70SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatDenseResetArray_C",NULL);CHKERRQ(ierr); 5958baccfbdSHong Zhang #if defined(PETSC_HAVE_ELEMENTAL) 5968baccfbdSHong Zhang ierr = PetscObjectComposeFunction((PetscObject)mat,"MatConvert_mpidense_elemental_C",NULL);CHKERRQ(ierr); 5978baccfbdSHong Zhang #endif 598bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMPIDenseSetPreallocation_C",NULL);CHKERRQ(ierr); 599bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMatMult_mpiaij_mpidense_C",NULL);CHKERRQ(ierr); 600bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMatMultSymbolic_mpiaij_mpidense_C",NULL);CHKERRQ(ierr); 601bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMatMultNumeric_mpiaij_mpidense_C",NULL);CHKERRQ(ierr); 60252c5f739Sprj- ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMatMult_nest_mpidense_C",NULL);CHKERRQ(ierr); 60352c5f739Sprj- ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMatMultSymbolic_nest_mpidense_C",NULL);CHKERRQ(ierr); 60452c5f739Sprj- ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMatMultNumeric_nest_mpidense_C",NULL);CHKERRQ(ierr); 6058baccfbdSHong Zhang ierr = PetscObjectComposeFunction((PetscObject)mat,"MatTransposeMatMult_mpiaij_mpidense_C",NULL);CHKERRQ(ierr); 6068baccfbdSHong Zhang ierr = PetscObjectComposeFunction((PetscObject)mat,"MatTransposeMatMultSymbolic_mpiaij_mpidense_C",NULL);CHKERRQ(ierr); 6078baccfbdSHong Zhang ierr = PetscObjectComposeFunction((PetscObject)mat,"MatTransposeMatMultNumeric_mpiaij_mpidense_C",NULL);CHKERRQ(ierr); 60886aefd0dSHong Zhang ierr = PetscObjectComposeFunction((PetscObject)mat,"MatDenseGetColumn_C",NULL);CHKERRQ(ierr); 60986aefd0dSHong Zhang ierr = PetscObjectComposeFunction((PetscObject)mat,"MatDenseRestoreColumn_C",NULL);CHKERRQ(ierr); 6103a40ed3dSBarry Smith PetscFunctionReturn(0); 6118965ea79SLois Curfman McInnes } 61239ddd567SLois Curfman McInnes 61352c5f739Sprj- PETSC_INTERN PetscErrorCode MatView_SeqDense(Mat,PetscViewer); 61452c5f739Sprj- 6159804daf3SBarry Smith #include <petscdraw.h> 6166849ba73SBarry Smith static PetscErrorCode MatView_MPIDense_ASCIIorDraworSocket(Mat mat,PetscViewer viewer) 6178965ea79SLois Curfman McInnes { 61839ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 619dfbe8321SBarry Smith PetscErrorCode ierr; 6207da1fb6eSBarry Smith PetscMPIInt rank = mdn->rank; 62119fd82e9SBarry Smith PetscViewerType vtype; 622ace3abfcSBarry Smith PetscBool iascii,isdraw; 623b0a32e0cSBarry Smith PetscViewer sviewer; 624f3ef73ceSBarry Smith PetscViewerFormat format; 6258965ea79SLois Curfman McInnes 6263a40ed3dSBarry Smith PetscFunctionBegin; 627251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 628251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr); 62932077d6dSBarry Smith if (iascii) { 630b0a32e0cSBarry Smith ierr = PetscViewerGetType(viewer,&vtype);CHKERRQ(ierr); 631b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 632456192e2SBarry Smith if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) { 6334e220ebcSLois Curfman McInnes MatInfo info; 634888f2ed8SSatish Balay ierr = MatGetInfo(mat,MAT_LOCAL,&info);CHKERRQ(ierr); 6351575c14dSBarry Smith ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr); 6367b23a99aSBarry 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); 637b0a32e0cSBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 6381575c14dSBarry Smith ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr); 6395d00a290SHong Zhang ierr = VecScatterView(mdn->Mvctx,viewer);CHKERRQ(ierr); 6403a40ed3dSBarry Smith PetscFunctionReturn(0); 641fb9695e5SSatish Balay } else if (format == PETSC_VIEWER_ASCII_INFO) { 6423a40ed3dSBarry Smith PetscFunctionReturn(0); 6438965ea79SLois Curfman McInnes } 644f1af5d2fSBarry Smith } else if (isdraw) { 645b0a32e0cSBarry Smith PetscDraw draw; 646ace3abfcSBarry Smith PetscBool isnull; 647f1af5d2fSBarry Smith 648b0a32e0cSBarry Smith ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 649b0a32e0cSBarry Smith ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); 650f1af5d2fSBarry Smith if (isnull) PetscFunctionReturn(0); 651f1af5d2fSBarry Smith } 65277ed5343SBarry Smith 6537da1fb6eSBarry Smith { 6548965ea79SLois Curfman McInnes /* assemble the entire matrix onto first processor. */ 6558965ea79SLois Curfman McInnes Mat A; 656d0f46423SBarry Smith PetscInt M = mat->rmap->N,N = mat->cmap->N,m,row,i,nz; 657ba8c8a56SBarry Smith PetscInt *cols; 658ba8c8a56SBarry Smith PetscScalar *vals; 6598965ea79SLois Curfman McInnes 660ce94432eSBarry Smith ierr = MatCreate(PetscObjectComm((PetscObject)mat),&A);CHKERRQ(ierr); 6618965ea79SLois Curfman McInnes if (!rank) { 662f69a0ea3SMatthew Knepley ierr = MatSetSizes(A,M,N,M,N);CHKERRQ(ierr); 6633a40ed3dSBarry Smith } else { 664f69a0ea3SMatthew Knepley ierr = MatSetSizes(A,0,0,M,N);CHKERRQ(ierr); 6658965ea79SLois Curfman McInnes } 6667adad957SLisandro Dalcin /* Since this is a temporary matrix, MATMPIDENSE instead of ((PetscObject)A)->type_name here is probably acceptable. */ 667878740d9SKris Buschelman ierr = MatSetType(A,MATMPIDENSE);CHKERRQ(ierr); 6680298fd71SBarry Smith ierr = MatMPIDenseSetPreallocation(A,NULL);CHKERRQ(ierr); 6693bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)A);CHKERRQ(ierr); 6708965ea79SLois Curfman McInnes 67139ddd567SLois Curfman McInnes /* Copy the matrix ... This isn't the most efficient means, 67239ddd567SLois Curfman McInnes but it's quick for now */ 67351022da4SBarry Smith A->insertmode = INSERT_VALUES; 6742205254eSKarl Rupp 6752205254eSKarl Rupp row = mat->rmap->rstart; 6762205254eSKarl Rupp m = mdn->A->rmap->n; 6778965ea79SLois Curfman McInnes for (i=0; i<m; i++) { 678ba8c8a56SBarry Smith ierr = MatGetRow_MPIDense(mat,row,&nz,&cols,&vals);CHKERRQ(ierr); 679ba8c8a56SBarry Smith ierr = MatSetValues_MPIDense(A,1,&row,nz,cols,vals,INSERT_VALUES);CHKERRQ(ierr); 680ba8c8a56SBarry Smith ierr = MatRestoreRow_MPIDense(mat,row,&nz,&cols,&vals);CHKERRQ(ierr); 68139ddd567SLois Curfman McInnes row++; 6828965ea79SLois Curfman McInnes } 6838965ea79SLois Curfman McInnes 6846d4a8577SBarry Smith ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 6856d4a8577SBarry Smith ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 6863f08860eSBarry Smith ierr = PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr); 687b9b97703SBarry Smith if (!rank) { 6881a9d3c3cSBarry Smith ierr = PetscObjectSetName((PetscObject)((Mat_MPIDense*)(A->data))->A,((PetscObject)mat)->name);CHKERRQ(ierr); 6897da1fb6eSBarry Smith ierr = MatView_SeqDense(((Mat_MPIDense*)(A->data))->A,sviewer);CHKERRQ(ierr); 6908965ea79SLois Curfman McInnes } 6913f08860eSBarry Smith ierr = PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr); 692b0a32e0cSBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 6936bf464f9SBarry Smith ierr = MatDestroy(&A);CHKERRQ(ierr); 6948965ea79SLois Curfman McInnes } 6953a40ed3dSBarry Smith PetscFunctionReturn(0); 6968965ea79SLois Curfman McInnes } 6978965ea79SLois Curfman McInnes 698dfbe8321SBarry Smith PetscErrorCode MatView_MPIDense(Mat mat,PetscViewer viewer) 6998965ea79SLois Curfman McInnes { 700dfbe8321SBarry Smith PetscErrorCode ierr; 701ace3abfcSBarry Smith PetscBool iascii,isbinary,isdraw,issocket; 7028965ea79SLois Curfman McInnes 703433994e6SBarry Smith PetscFunctionBegin; 704251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 705251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr); 706251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSOCKET,&issocket);CHKERRQ(ierr); 707251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr); 7080f5bd95cSBarry Smith 70932077d6dSBarry Smith if (iascii || issocket || isdraw) { 710f1af5d2fSBarry Smith ierr = MatView_MPIDense_ASCIIorDraworSocket(mat,viewer);CHKERRQ(ierr); 7110f5bd95cSBarry Smith } else if (isbinary) { 712*8491ab44SLisandro Dalcin ierr = MatView_Dense_Binary(mat,viewer);CHKERRQ(ierr); 71311aeaf0aSBarry Smith } 7143a40ed3dSBarry Smith PetscFunctionReturn(0); 7158965ea79SLois Curfman McInnes } 7168965ea79SLois Curfman McInnes 717dfbe8321SBarry Smith PetscErrorCode MatGetInfo_MPIDense(Mat A,MatInfoType flag,MatInfo *info) 7188965ea79SLois Curfman McInnes { 7193501a2bdSLois Curfman McInnes Mat_MPIDense *mat = (Mat_MPIDense*)A->data; 7203501a2bdSLois Curfman McInnes Mat mdn = mat->A; 721dfbe8321SBarry Smith PetscErrorCode ierr; 7223966268fSBarry Smith PetscLogDouble isend[5],irecv[5]; 7238965ea79SLois Curfman McInnes 7243a40ed3dSBarry Smith PetscFunctionBegin; 7254e220ebcSLois Curfman McInnes info->block_size = 1.0; 7262205254eSKarl Rupp 7274e220ebcSLois Curfman McInnes ierr = MatGetInfo(mdn,MAT_LOCAL,info);CHKERRQ(ierr); 7282205254eSKarl Rupp 7294e220ebcSLois Curfman McInnes isend[0] = info->nz_used; isend[1] = info->nz_allocated; isend[2] = info->nz_unneeded; 7304e220ebcSLois Curfman McInnes isend[3] = info->memory; isend[4] = info->mallocs; 7318965ea79SLois Curfman McInnes if (flag == MAT_LOCAL) { 7324e220ebcSLois Curfman McInnes info->nz_used = isend[0]; 7334e220ebcSLois Curfman McInnes info->nz_allocated = isend[1]; 7344e220ebcSLois Curfman McInnes info->nz_unneeded = isend[2]; 7354e220ebcSLois Curfman McInnes info->memory = isend[3]; 7364e220ebcSLois Curfman McInnes info->mallocs = isend[4]; 7378965ea79SLois Curfman McInnes } else if (flag == MAT_GLOBAL_MAX) { 7383966268fSBarry Smith ierr = MPIU_Allreduce(isend,irecv,5,MPIU_PETSCLOGDOUBLE,MPI_MAX,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 7392205254eSKarl Rupp 7404e220ebcSLois Curfman McInnes info->nz_used = irecv[0]; 7414e220ebcSLois Curfman McInnes info->nz_allocated = irecv[1]; 7424e220ebcSLois Curfman McInnes info->nz_unneeded = irecv[2]; 7434e220ebcSLois Curfman McInnes info->memory = irecv[3]; 7444e220ebcSLois Curfman McInnes info->mallocs = irecv[4]; 7458965ea79SLois Curfman McInnes } else if (flag == MAT_GLOBAL_SUM) { 7463966268fSBarry Smith ierr = MPIU_Allreduce(isend,irecv,5,MPIU_PETSCLOGDOUBLE,MPI_SUM,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 7472205254eSKarl Rupp 7484e220ebcSLois Curfman McInnes info->nz_used = irecv[0]; 7494e220ebcSLois Curfman McInnes info->nz_allocated = irecv[1]; 7504e220ebcSLois Curfman McInnes info->nz_unneeded = irecv[2]; 7514e220ebcSLois Curfman McInnes info->memory = irecv[3]; 7524e220ebcSLois Curfman McInnes info->mallocs = irecv[4]; 7538965ea79SLois Curfman McInnes } 7544e220ebcSLois Curfman McInnes info->fill_ratio_given = 0; /* no parallel LU/ILU/Cholesky */ 7554e220ebcSLois Curfman McInnes info->fill_ratio_needed = 0; 7564e220ebcSLois Curfman McInnes info->factor_mallocs = 0; 7573a40ed3dSBarry Smith PetscFunctionReturn(0); 7588965ea79SLois Curfman McInnes } 7598965ea79SLois Curfman McInnes 760ace3abfcSBarry Smith PetscErrorCode MatSetOption_MPIDense(Mat A,MatOption op,PetscBool flg) 7618965ea79SLois Curfman McInnes { 76239ddd567SLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense*)A->data; 763dfbe8321SBarry Smith PetscErrorCode ierr; 7648965ea79SLois Curfman McInnes 7653a40ed3dSBarry Smith PetscFunctionBegin; 76612c028f9SKris Buschelman switch (op) { 767512a5fc5SBarry Smith case MAT_NEW_NONZERO_LOCATIONS: 76812c028f9SKris Buschelman case MAT_NEW_NONZERO_LOCATION_ERR: 76912c028f9SKris Buschelman case MAT_NEW_NONZERO_ALLOCATION_ERR: 77043674050SBarry Smith MatCheckPreallocated(A,1); 7714e0d8c25SBarry Smith ierr = MatSetOption(a->A,op,flg);CHKERRQ(ierr); 77212c028f9SKris Buschelman break; 77312c028f9SKris Buschelman case MAT_ROW_ORIENTED: 77443674050SBarry Smith MatCheckPreallocated(A,1); 7754e0d8c25SBarry Smith a->roworiented = flg; 7764e0d8c25SBarry Smith ierr = MatSetOption(a->A,op,flg);CHKERRQ(ierr); 77712c028f9SKris Buschelman break; 7784e0d8c25SBarry Smith case MAT_NEW_DIAGONALS: 77913fa8e87SLisandro Dalcin case MAT_KEEP_NONZERO_PATTERN: 78012c028f9SKris Buschelman case MAT_USE_HASH_TABLE: 781071fcb05SBarry Smith case MAT_SORTED_FULL: 782290bbb0aSBarry Smith ierr = PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);CHKERRQ(ierr); 78312c028f9SKris Buschelman break; 78412c028f9SKris Buschelman case MAT_IGNORE_OFF_PROC_ENTRIES: 7854e0d8c25SBarry Smith a->donotstash = flg; 78612c028f9SKris Buschelman break; 78777e54ba9SKris Buschelman case MAT_SYMMETRIC: 78877e54ba9SKris Buschelman case MAT_STRUCTURALLY_SYMMETRIC: 7899a4540c5SBarry Smith case MAT_HERMITIAN: 7909a4540c5SBarry Smith case MAT_SYMMETRY_ETERNAL: 791600fe468SBarry Smith case MAT_IGNORE_LOWER_TRIANGULAR: 7925d7aebe8SStefano Zampini case MAT_IGNORE_ZERO_ENTRIES: 793290bbb0aSBarry Smith ierr = PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);CHKERRQ(ierr); 79477e54ba9SKris Buschelman break; 79512c028f9SKris Buschelman default: 796e32f2f54SBarry Smith SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"unknown option %s",MatOptions[op]); 7973a40ed3dSBarry Smith } 7983a40ed3dSBarry Smith PetscFunctionReturn(0); 7998965ea79SLois Curfman McInnes } 8008965ea79SLois Curfman McInnes 8018965ea79SLois Curfman McInnes 802dfbe8321SBarry Smith PetscErrorCode MatDiagonalScale_MPIDense(Mat A,Vec ll,Vec rr) 8035b2fa520SLois Curfman McInnes { 8045b2fa520SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)A->data; 8055b2fa520SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)mdn->A->data; 806bca11509SBarry Smith const PetscScalar *l,*r; 807bca11509SBarry Smith PetscScalar x,*v; 808dfbe8321SBarry Smith PetscErrorCode ierr; 809d0f46423SBarry Smith PetscInt i,j,s2a,s3a,s2,s3,m=mdn->A->rmap->n,n=mdn->A->cmap->n; 8105b2fa520SLois Curfman McInnes 8115b2fa520SLois Curfman McInnes PetscFunctionBegin; 81272d926a5SLois Curfman McInnes ierr = MatGetLocalSize(A,&s2,&s3);CHKERRQ(ierr); 8135b2fa520SLois Curfman McInnes if (ll) { 81472d926a5SLois Curfman McInnes ierr = VecGetLocalSize(ll,&s2a);CHKERRQ(ierr); 815e32f2f54SBarry Smith if (s2a != s2) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Left scaling vector non-conforming local size, %d != %d.", s2a, s2); 816bca11509SBarry Smith ierr = VecGetArrayRead(ll,&l);CHKERRQ(ierr); 8175b2fa520SLois Curfman McInnes for (i=0; i<m; i++) { 8185b2fa520SLois Curfman McInnes x = l[i]; 8195b2fa520SLois Curfman McInnes v = mat->v + i; 8205b2fa520SLois Curfman McInnes for (j=0; j<n; j++) { (*v) *= x; v+= m;} 8215b2fa520SLois Curfman McInnes } 822bca11509SBarry Smith ierr = VecRestoreArrayRead(ll,&l);CHKERRQ(ierr); 823efee365bSSatish Balay ierr = PetscLogFlops(n*m);CHKERRQ(ierr); 8245b2fa520SLois Curfman McInnes } 8255b2fa520SLois Curfman McInnes if (rr) { 826175be7b4SMatthew Knepley ierr = VecGetLocalSize(rr,&s3a);CHKERRQ(ierr); 827e32f2f54SBarry Smith if (s3a != s3) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Right scaling vec non-conforming local size, %d != %d.", s3a, s3); 828ca9f406cSSatish Balay ierr = VecScatterBegin(mdn->Mvctx,rr,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 829ca9f406cSSatish Balay ierr = VecScatterEnd(mdn->Mvctx,rr,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 830bca11509SBarry Smith ierr = VecGetArrayRead(mdn->lvec,&r);CHKERRQ(ierr); 8315b2fa520SLois Curfman McInnes for (i=0; i<n; i++) { 8325b2fa520SLois Curfman McInnes x = r[i]; 8335b2fa520SLois Curfman McInnes v = mat->v + i*m; 8342205254eSKarl Rupp for (j=0; j<m; j++) (*v++) *= x; 8355b2fa520SLois Curfman McInnes } 836bca11509SBarry Smith ierr = VecRestoreArrayRead(mdn->lvec,&r);CHKERRQ(ierr); 837efee365bSSatish Balay ierr = PetscLogFlops(n*m);CHKERRQ(ierr); 8385b2fa520SLois Curfman McInnes } 8395b2fa520SLois Curfman McInnes PetscFunctionReturn(0); 8405b2fa520SLois Curfman McInnes } 8415b2fa520SLois Curfman McInnes 842dfbe8321SBarry Smith PetscErrorCode MatNorm_MPIDense(Mat A,NormType type,PetscReal *nrm) 843096963f5SLois Curfman McInnes { 8443501a2bdSLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense*)A->data; 8453501a2bdSLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)mdn->A->data; 846dfbe8321SBarry Smith PetscErrorCode ierr; 84713f74950SBarry Smith PetscInt i,j; 848329f5518SBarry Smith PetscReal sum = 0.0; 84987828ca2SBarry Smith PetscScalar *v = mat->v; 8503501a2bdSLois Curfman McInnes 8513a40ed3dSBarry Smith PetscFunctionBegin; 8523501a2bdSLois Curfman McInnes if (mdn->size == 1) { 853064f8208SBarry Smith ierr = MatNorm(mdn->A,type,nrm);CHKERRQ(ierr); 8543501a2bdSLois Curfman McInnes } else { 8553501a2bdSLois Curfman McInnes if (type == NORM_FROBENIUS) { 856d0f46423SBarry Smith for (i=0; i<mdn->A->cmap->n*mdn->A->rmap->n; i++) { 857329f5518SBarry Smith sum += PetscRealPart(PetscConj(*v)*(*v)); v++; 8583501a2bdSLois Curfman McInnes } 859b2566f29SBarry Smith ierr = MPIU_Allreduce(&sum,nrm,1,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 8608f1a2a5eSBarry Smith *nrm = PetscSqrtReal(*nrm); 861dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*mdn->A->cmap->n*mdn->A->rmap->n);CHKERRQ(ierr); 8623a40ed3dSBarry Smith } else if (type == NORM_1) { 863329f5518SBarry Smith PetscReal *tmp,*tmp2; 864580bdb30SBarry Smith ierr = PetscCalloc2(A->cmap->N,&tmp,A->cmap->N,&tmp2);CHKERRQ(ierr); 865064f8208SBarry Smith *nrm = 0.0; 8663501a2bdSLois Curfman McInnes v = mat->v; 867d0f46423SBarry Smith for (j=0; j<mdn->A->cmap->n; j++) { 868d0f46423SBarry Smith for (i=0; i<mdn->A->rmap->n; i++) { 86967e560aaSBarry Smith tmp[j] += PetscAbsScalar(*v); v++; 8703501a2bdSLois Curfman McInnes } 8713501a2bdSLois Curfman McInnes } 872b2566f29SBarry Smith ierr = MPIU_Allreduce(tmp,tmp2,A->cmap->N,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 873d0f46423SBarry Smith for (j=0; j<A->cmap->N; j++) { 874064f8208SBarry Smith if (tmp2[j] > *nrm) *nrm = tmp2[j]; 8753501a2bdSLois Curfman McInnes } 8768627564fSBarry Smith ierr = PetscFree2(tmp,tmp2);CHKERRQ(ierr); 877d0f46423SBarry Smith ierr = PetscLogFlops(A->cmap->n*A->rmap->n);CHKERRQ(ierr); 8783a40ed3dSBarry Smith } else if (type == NORM_INFINITY) { /* max row norm */ 879329f5518SBarry Smith PetscReal ntemp; 8803501a2bdSLois Curfman McInnes ierr = MatNorm(mdn->A,type,&ntemp);CHKERRQ(ierr); 881b2566f29SBarry Smith ierr = MPIU_Allreduce(&ntemp,nrm,1,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 882ce94432eSBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support for two norm"); 8833501a2bdSLois Curfman McInnes } 8843a40ed3dSBarry Smith PetscFunctionReturn(0); 8853501a2bdSLois Curfman McInnes } 8863501a2bdSLois Curfman McInnes 887fc4dec0aSBarry Smith PetscErrorCode MatTranspose_MPIDense(Mat A,MatReuse reuse,Mat *matout) 8883501a2bdSLois Curfman McInnes { 8893501a2bdSLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense*)A->data; 8903501a2bdSLois Curfman McInnes Mat_SeqDense *Aloc = (Mat_SeqDense*)a->A->data; 8913501a2bdSLois Curfman McInnes Mat B; 892d0f46423SBarry Smith PetscInt M = A->rmap->N,N = A->cmap->N,m,n,*rwork,rstart = A->rmap->rstart; 8936849ba73SBarry Smith PetscErrorCode ierr; 89413f74950SBarry Smith PetscInt j,i; 89587828ca2SBarry Smith PetscScalar *v; 8963501a2bdSLois Curfman McInnes 8973a40ed3dSBarry Smith PetscFunctionBegin; 898cf37664fSBarry Smith if (reuse == MAT_INITIAL_MATRIX || reuse == MAT_INPLACE_MATRIX) { 899ce94432eSBarry Smith ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr); 900d0f46423SBarry Smith ierr = MatSetSizes(B,A->cmap->n,A->rmap->n,N,M);CHKERRQ(ierr); 9017adad957SLisandro Dalcin ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr); 9020298fd71SBarry Smith ierr = MatMPIDenseSetPreallocation(B,NULL);CHKERRQ(ierr); 903fc4dec0aSBarry Smith } else { 904fc4dec0aSBarry Smith B = *matout; 905fc4dec0aSBarry Smith } 9063501a2bdSLois Curfman McInnes 907d0f46423SBarry Smith m = a->A->rmap->n; n = a->A->cmap->n; v = Aloc->v; 908785e854fSJed Brown ierr = PetscMalloc1(m,&rwork);CHKERRQ(ierr); 9093501a2bdSLois Curfman McInnes for (i=0; i<m; i++) rwork[i] = rstart + i; 9101acff37aSSatish Balay for (j=0; j<n; j++) { 9113501a2bdSLois Curfman McInnes ierr = MatSetValues(B,1,&j,m,rwork,v,INSERT_VALUES);CHKERRQ(ierr); 9123501a2bdSLois Curfman McInnes v += m; 9133501a2bdSLois Curfman McInnes } 914606d414cSSatish Balay ierr = PetscFree(rwork);CHKERRQ(ierr); 9156d4a8577SBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 9166d4a8577SBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 917cf37664fSBarry Smith if (reuse == MAT_INITIAL_MATRIX || reuse == MAT_REUSE_MATRIX) { 9183501a2bdSLois Curfman McInnes *matout = B; 9193501a2bdSLois Curfman McInnes } else { 92028be2f97SBarry Smith ierr = MatHeaderMerge(A,&B);CHKERRQ(ierr); 9213501a2bdSLois Curfman McInnes } 9223a40ed3dSBarry Smith PetscFunctionReturn(0); 923096963f5SLois Curfman McInnes } 924096963f5SLois Curfman McInnes 9256849ba73SBarry Smith static PetscErrorCode MatDuplicate_MPIDense(Mat,MatDuplicateOption,Mat*); 92652c5f739Sprj- PETSC_INTERN PetscErrorCode MatScale_MPIDense(Mat,PetscScalar); 9278965ea79SLois Curfman McInnes 9284994cf47SJed Brown PetscErrorCode MatSetUp_MPIDense(Mat A) 929273d9f13SBarry Smith { 930dfbe8321SBarry Smith PetscErrorCode ierr; 931273d9f13SBarry Smith 932273d9f13SBarry Smith PetscFunctionBegin; 933273d9f13SBarry Smith ierr = MatMPIDenseSetPreallocation(A,0);CHKERRQ(ierr); 934273d9f13SBarry Smith PetscFunctionReturn(0); 935273d9f13SBarry Smith } 936273d9f13SBarry Smith 937488007eeSBarry Smith PetscErrorCode MatAXPY_MPIDense(Mat Y,PetscScalar alpha,Mat X,MatStructure str) 938488007eeSBarry Smith { 939488007eeSBarry Smith PetscErrorCode ierr; 940488007eeSBarry Smith Mat_MPIDense *A = (Mat_MPIDense*)Y->data, *B = (Mat_MPIDense*)X->data; 941488007eeSBarry Smith 942488007eeSBarry Smith PetscFunctionBegin; 943488007eeSBarry Smith ierr = MatAXPY(A->A,alpha,B->A,str);CHKERRQ(ierr); 944a3fa217bSJose E. Roman ierr = PetscObjectStateIncrease((PetscObject)Y);CHKERRQ(ierr); 945488007eeSBarry Smith PetscFunctionReturn(0); 946488007eeSBarry Smith } 947488007eeSBarry Smith 9487087cfbeSBarry Smith PetscErrorCode MatConjugate_MPIDense(Mat mat) 949ba337c44SJed Brown { 950ba337c44SJed Brown Mat_MPIDense *a = (Mat_MPIDense*)mat->data; 951ba337c44SJed Brown PetscErrorCode ierr; 952ba337c44SJed Brown 953ba337c44SJed Brown PetscFunctionBegin; 954ba337c44SJed Brown ierr = MatConjugate(a->A);CHKERRQ(ierr); 955ba337c44SJed Brown PetscFunctionReturn(0); 956ba337c44SJed Brown } 957ba337c44SJed Brown 958ba337c44SJed Brown PetscErrorCode MatRealPart_MPIDense(Mat A) 959ba337c44SJed Brown { 960ba337c44SJed Brown Mat_MPIDense *a = (Mat_MPIDense*)A->data; 961ba337c44SJed Brown PetscErrorCode ierr; 962ba337c44SJed Brown 963ba337c44SJed Brown PetscFunctionBegin; 964ba337c44SJed Brown ierr = MatRealPart(a->A);CHKERRQ(ierr); 965ba337c44SJed Brown PetscFunctionReturn(0); 966ba337c44SJed Brown } 967ba337c44SJed Brown 968ba337c44SJed Brown PetscErrorCode MatImaginaryPart_MPIDense(Mat A) 969ba337c44SJed Brown { 970ba337c44SJed Brown Mat_MPIDense *a = (Mat_MPIDense*)A->data; 971ba337c44SJed Brown PetscErrorCode ierr; 972ba337c44SJed Brown 973ba337c44SJed Brown PetscFunctionBegin; 974ba337c44SJed Brown ierr = MatImaginaryPart(a->A);CHKERRQ(ierr); 975ba337c44SJed Brown PetscFunctionReturn(0); 976ba337c44SJed Brown } 977ba337c44SJed Brown 97849a6ff4bSBarry Smith static PetscErrorCode MatGetColumnVector_MPIDense(Mat A,Vec v,PetscInt col) 97949a6ff4bSBarry Smith { 98049a6ff4bSBarry Smith PetscErrorCode ierr; 98149a6ff4bSBarry Smith PetscScalar *x; 98249a6ff4bSBarry Smith const PetscScalar *a; 98349a6ff4bSBarry Smith PetscInt lda; 98449a6ff4bSBarry Smith 98549a6ff4bSBarry Smith PetscFunctionBegin; 98649a6ff4bSBarry Smith if (A->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 98749a6ff4bSBarry Smith ierr = MatDenseGetLDA(A,&lda);CHKERRQ(ierr); 98849a6ff4bSBarry Smith ierr = MatDenseGetArrayRead(A,&a);CHKERRQ(ierr); 98949a6ff4bSBarry Smith ierr = VecGetArray(v,&x);CHKERRQ(ierr); 990580bdb30SBarry Smith ierr = PetscArraycpy(x,a+col*lda,A->rmap->n);CHKERRQ(ierr); 99149a6ff4bSBarry Smith ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 99249a6ff4bSBarry Smith ierr = MatDenseGetArrayRead(A,&a);CHKERRQ(ierr); 99349a6ff4bSBarry Smith PetscFunctionReturn(0); 99449a6ff4bSBarry Smith } 99549a6ff4bSBarry Smith 99652c5f739Sprj- PETSC_INTERN PetscErrorCode MatGetColumnNorms_SeqDense(Mat,NormType,PetscReal*); 99752c5f739Sprj- 9980716a85fSBarry Smith PetscErrorCode MatGetColumnNorms_MPIDense(Mat A,NormType type,PetscReal *norms) 9990716a85fSBarry Smith { 10000716a85fSBarry Smith PetscErrorCode ierr; 10010716a85fSBarry Smith PetscInt i,n; 10020716a85fSBarry Smith Mat_MPIDense *a = (Mat_MPIDense*) A->data; 10030716a85fSBarry Smith PetscReal *work; 10040716a85fSBarry Smith 10050716a85fSBarry Smith PetscFunctionBegin; 10060298fd71SBarry Smith ierr = MatGetSize(A,NULL,&n);CHKERRQ(ierr); 1007785e854fSJed Brown ierr = PetscMalloc1(n,&work);CHKERRQ(ierr); 10080716a85fSBarry Smith ierr = MatGetColumnNorms_SeqDense(a->A,type,work);CHKERRQ(ierr); 10090716a85fSBarry Smith if (type == NORM_2) { 10100716a85fSBarry Smith for (i=0; i<n; i++) work[i] *= work[i]; 10110716a85fSBarry Smith } 10120716a85fSBarry Smith if (type == NORM_INFINITY) { 1013b2566f29SBarry Smith ierr = MPIU_Allreduce(work,norms,n,MPIU_REAL,MPIU_MAX,A->hdr.comm);CHKERRQ(ierr); 10140716a85fSBarry Smith } else { 1015b2566f29SBarry Smith ierr = MPIU_Allreduce(work,norms,n,MPIU_REAL,MPIU_SUM,A->hdr.comm);CHKERRQ(ierr); 10160716a85fSBarry Smith } 10170716a85fSBarry Smith ierr = PetscFree(work);CHKERRQ(ierr); 10180716a85fSBarry Smith if (type == NORM_2) { 10198f1a2a5eSBarry Smith for (i=0; i<n; i++) norms[i] = PetscSqrtReal(norms[i]); 10200716a85fSBarry Smith } 10210716a85fSBarry Smith PetscFunctionReturn(0); 10220716a85fSBarry Smith } 10230716a85fSBarry Smith 102473a71a0fSBarry Smith static PetscErrorCode MatSetRandom_MPIDense(Mat x,PetscRandom rctx) 102573a71a0fSBarry Smith { 102673a71a0fSBarry Smith Mat_MPIDense *d = (Mat_MPIDense*)x->data; 102773a71a0fSBarry Smith PetscErrorCode ierr; 102873a71a0fSBarry Smith PetscScalar *a; 102973a71a0fSBarry Smith PetscInt m,n,i; 103073a71a0fSBarry Smith 103173a71a0fSBarry Smith PetscFunctionBegin; 103273a71a0fSBarry Smith ierr = MatGetSize(d->A,&m,&n);CHKERRQ(ierr); 10338c778c55SBarry Smith ierr = MatDenseGetArray(d->A,&a);CHKERRQ(ierr); 103473a71a0fSBarry Smith for (i=0; i<m*n; i++) { 103573a71a0fSBarry Smith ierr = PetscRandomGetValue(rctx,a+i);CHKERRQ(ierr); 103673a71a0fSBarry Smith } 10378c778c55SBarry Smith ierr = MatDenseRestoreArray(d->A,&a);CHKERRQ(ierr); 103873a71a0fSBarry Smith PetscFunctionReturn(0); 103973a71a0fSBarry Smith } 104073a71a0fSBarry Smith 104152c5f739Sprj- PETSC_INTERN PetscErrorCode MatMatMultNumeric_MPIDense(Mat A,Mat,Mat); 1042fd4e9aacSBarry Smith 10433b49f96aSBarry Smith static PetscErrorCode MatMissingDiagonal_MPIDense(Mat A,PetscBool *missing,PetscInt *d) 10443b49f96aSBarry Smith { 10453b49f96aSBarry Smith PetscFunctionBegin; 10463b49f96aSBarry Smith *missing = PETSC_FALSE; 10473b49f96aSBarry Smith PetscFunctionReturn(0); 10483b49f96aSBarry Smith } 10493b49f96aSBarry Smith 1050cc48ffa7SToby Isaac static PetscErrorCode MatMatTransposeMult_MPIDense_MPIDense(Mat,Mat,MatReuse,PetscReal,Mat*); 1051cc48ffa7SToby Isaac static PetscErrorCode MatMatTransposeMultSymbolic_MPIDense_MPIDense(Mat,Mat,PetscReal,Mat*); 1052cc48ffa7SToby Isaac static PetscErrorCode MatMatTransposeMultNumeric_MPIDense_MPIDense(Mat,Mat,Mat); 1053cc48ffa7SToby Isaac 10548965ea79SLois Curfman McInnes /* -------------------------------------------------------------------*/ 105509dc0095SBarry Smith static struct _MatOps MatOps_Values = { MatSetValues_MPIDense, 105609dc0095SBarry Smith MatGetRow_MPIDense, 105709dc0095SBarry Smith MatRestoreRow_MPIDense, 105809dc0095SBarry Smith MatMult_MPIDense, 105997304618SKris Buschelman /* 4*/ MatMultAdd_MPIDense, 10607c922b88SBarry Smith MatMultTranspose_MPIDense, 10617c922b88SBarry Smith MatMultTransposeAdd_MPIDense, 10628965ea79SLois Curfman McInnes 0, 106309dc0095SBarry Smith 0, 106409dc0095SBarry Smith 0, 106597304618SKris Buschelman /* 10*/ 0, 106609dc0095SBarry Smith 0, 106709dc0095SBarry Smith 0, 106809dc0095SBarry Smith 0, 106909dc0095SBarry Smith MatTranspose_MPIDense, 107097304618SKris Buschelman /* 15*/ MatGetInfo_MPIDense, 10716e4ee0c6SHong Zhang MatEqual_MPIDense, 107209dc0095SBarry Smith MatGetDiagonal_MPIDense, 10735b2fa520SLois Curfman McInnes MatDiagonalScale_MPIDense, 107409dc0095SBarry Smith MatNorm_MPIDense, 107597304618SKris Buschelman /* 20*/ MatAssemblyBegin_MPIDense, 107609dc0095SBarry Smith MatAssemblyEnd_MPIDense, 107709dc0095SBarry Smith MatSetOption_MPIDense, 107809dc0095SBarry Smith MatZeroEntries_MPIDense, 1079d519adbfSMatthew Knepley /* 24*/ MatZeroRows_MPIDense, 1080919b68f7SBarry Smith 0, 108101b82886SBarry Smith 0, 108201b82886SBarry Smith 0, 108301b82886SBarry Smith 0, 10844994cf47SJed Brown /* 29*/ MatSetUp_MPIDense, 1085273d9f13SBarry Smith 0, 108609dc0095SBarry Smith 0, 1087c56a70eeSBarry Smith MatGetDiagonalBlock_MPIDense, 10888c778c55SBarry Smith 0, 1089d519adbfSMatthew Knepley /* 34*/ MatDuplicate_MPIDense, 109009dc0095SBarry Smith 0, 109109dc0095SBarry Smith 0, 109209dc0095SBarry Smith 0, 109309dc0095SBarry Smith 0, 1094d519adbfSMatthew Knepley /* 39*/ MatAXPY_MPIDense, 10957dae84e0SHong Zhang MatCreateSubMatrices_MPIDense, 109609dc0095SBarry Smith 0, 109709dc0095SBarry Smith MatGetValues_MPIDense, 109809dc0095SBarry Smith 0, 1099d519adbfSMatthew Knepley /* 44*/ 0, 110009dc0095SBarry Smith MatScale_MPIDense, 11017d68702bSBarry Smith MatShift_Basic, 110209dc0095SBarry Smith 0, 110309dc0095SBarry Smith 0, 110473a71a0fSBarry Smith /* 49*/ MatSetRandom_MPIDense, 110509dc0095SBarry Smith 0, 110609dc0095SBarry Smith 0, 110709dc0095SBarry Smith 0, 110809dc0095SBarry Smith 0, 1109d519adbfSMatthew Knepley /* 54*/ 0, 111009dc0095SBarry Smith 0, 111109dc0095SBarry Smith 0, 111209dc0095SBarry Smith 0, 111309dc0095SBarry Smith 0, 11147dae84e0SHong Zhang /* 59*/ MatCreateSubMatrix_MPIDense, 1115b9b97703SBarry Smith MatDestroy_MPIDense, 1116b9b97703SBarry Smith MatView_MPIDense, 1117357abbc8SBarry Smith 0, 111897304618SKris Buschelman 0, 1119d519adbfSMatthew Knepley /* 64*/ 0, 112097304618SKris Buschelman 0, 112197304618SKris Buschelman 0, 112297304618SKris Buschelman 0, 112397304618SKris Buschelman 0, 1124d519adbfSMatthew Knepley /* 69*/ 0, 112597304618SKris Buschelman 0, 112697304618SKris Buschelman 0, 112797304618SKris Buschelman 0, 112897304618SKris Buschelman 0, 1129d519adbfSMatthew Knepley /* 74*/ 0, 113097304618SKris Buschelman 0, 113197304618SKris Buschelman 0, 113297304618SKris Buschelman 0, 113397304618SKris Buschelman 0, 1134d519adbfSMatthew Knepley /* 79*/ 0, 113597304618SKris Buschelman 0, 113697304618SKris Buschelman 0, 113797304618SKris Buschelman 0, 11385bba2384SShri Abhyankar /* 83*/ MatLoad_MPIDense, 1139865e5f61SKris Buschelman 0, 1140865e5f61SKris Buschelman 0, 1141865e5f61SKris Buschelman 0, 1142865e5f61SKris Buschelman 0, 1143865e5f61SKris Buschelman 0, 11449775878aSHong Zhang #if defined(PETSC_HAVE_ELEMENTAL) 1145320f2790SHong Zhang /* 89*/ MatMatMult_MPIDense_MPIDense, 1146320f2790SHong Zhang MatMatMultSymbolic_MPIDense_MPIDense, 11479775878aSHong Zhang #else 11489775878aSHong Zhang /* 89*/ 0, 1149865e5f61SKris Buschelman 0, 11509775878aSHong Zhang #endif 1151fd4e9aacSBarry Smith MatMatMultNumeric_MPIDense, 11522fbe02b9SBarry Smith 0, 1153ba337c44SJed Brown 0, 1154d519adbfSMatthew Knepley /* 94*/ 0, 1155cc48ffa7SToby Isaac MatMatTransposeMult_MPIDense_MPIDense, 1156cc48ffa7SToby Isaac MatMatTransposeMultSymbolic_MPIDense_MPIDense, 1157cc48ffa7SToby Isaac MatMatTransposeMultNumeric_MPIDense_MPIDense, 1158ba337c44SJed Brown 0, 1159ba337c44SJed Brown /* 99*/ 0, 1160ba337c44SJed Brown 0, 1161ba337c44SJed Brown 0, 1162ba337c44SJed Brown MatConjugate_MPIDense, 1163ba337c44SJed Brown 0, 1164ba337c44SJed Brown /*104*/ 0, 1165ba337c44SJed Brown MatRealPart_MPIDense, 1166ba337c44SJed Brown MatImaginaryPart_MPIDense, 116786d161a7SShri Abhyankar 0, 116886d161a7SShri Abhyankar 0, 116986d161a7SShri Abhyankar /*109*/ 0, 117086d161a7SShri Abhyankar 0, 117186d161a7SShri Abhyankar 0, 117249a6ff4bSBarry Smith MatGetColumnVector_MPIDense, 11733b49f96aSBarry Smith MatMissingDiagonal_MPIDense, 117486d161a7SShri Abhyankar /*114*/ 0, 117586d161a7SShri Abhyankar 0, 117686d161a7SShri Abhyankar 0, 117786d161a7SShri Abhyankar 0, 117886d161a7SShri Abhyankar 0, 117986d161a7SShri Abhyankar /*119*/ 0, 118086d161a7SShri Abhyankar 0, 118186d161a7SShri Abhyankar 0, 11820716a85fSBarry Smith 0, 11830716a85fSBarry Smith 0, 11840716a85fSBarry Smith /*124*/ 0, 11853964eb88SJed Brown MatGetColumnNorms_MPIDense, 11863964eb88SJed Brown 0, 11873964eb88SJed Brown 0, 11883964eb88SJed Brown 0, 11893964eb88SJed Brown /*129*/ 0, 1190cb20be35SHong Zhang MatTransposeMatMult_MPIDense_MPIDense, 1191cb20be35SHong Zhang MatTransposeMatMultSymbolic_MPIDense_MPIDense, 1192cb20be35SHong Zhang MatTransposeMatMultNumeric_MPIDense_MPIDense, 11933964eb88SJed Brown 0, 11943964eb88SJed Brown /*134*/ 0, 11953964eb88SJed Brown 0, 11963964eb88SJed Brown 0, 11973964eb88SJed Brown 0, 11983964eb88SJed Brown 0, 11993964eb88SJed Brown /*139*/ 0, 1200f9426fe0SMark Adams 0, 120194e2cb23SJakub Kruzik 0, 120294e2cb23SJakub Kruzik 0, 120394e2cb23SJakub Kruzik 0, 12041ae5eee8SJakub Kruzik /*144*/ MatCreateMPIMatConcatenateSeqMat_MPIDense 1205ba337c44SJed Brown }; 12068965ea79SLois Curfman McInnes 12077087cfbeSBarry Smith PetscErrorCode MatMPIDenseSetPreallocation_MPIDense(Mat mat,PetscScalar *data) 1208a23d5eceSKris Buschelman { 1209a23d5eceSKris Buschelman Mat_MPIDense *a; 1210dfbe8321SBarry Smith PetscErrorCode ierr; 1211a23d5eceSKris Buschelman 1212a23d5eceSKris Buschelman PetscFunctionBegin; 1213430ada49SBarry Smith if (mat->preallocated) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Matrix has already been preallocated"); 1214a23d5eceSKris Buschelman mat->preallocated = PETSC_TRUE; 1215a23d5eceSKris Buschelman /* Note: For now, when data is specified above, this assumes the user correctly 1216a23d5eceSKris Buschelman allocates the local dense storage space. We should add error checking. */ 1217a23d5eceSKris Buschelman 1218a23d5eceSKris Buschelman a = (Mat_MPIDense*)mat->data; 121934ef9618SShri Abhyankar ierr = PetscLayoutSetUp(mat->rmap);CHKERRQ(ierr); 122034ef9618SShri Abhyankar ierr = PetscLayoutSetUp(mat->cmap);CHKERRQ(ierr); 122134ef9618SShri Abhyankar a->nvec = mat->cmap->n; 122234ef9618SShri Abhyankar 1223f69a0ea3SMatthew Knepley ierr = MatCreate(PETSC_COMM_SELF,&a->A);CHKERRQ(ierr); 1224d0f46423SBarry Smith ierr = MatSetSizes(a->A,mat->rmap->n,mat->cmap->N,mat->rmap->n,mat->cmap->N);CHKERRQ(ierr); 12255c5985e7SKris Buschelman ierr = MatSetType(a->A,MATSEQDENSE);CHKERRQ(ierr); 12265c5985e7SKris Buschelman ierr = MatSeqDenseSetPreallocation(a->A,data);CHKERRQ(ierr); 12273bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)a->A);CHKERRQ(ierr); 1228a23d5eceSKris Buschelman PetscFunctionReturn(0); 1229a23d5eceSKris Buschelman } 1230a23d5eceSKris Buschelman 123165b80a83SHong Zhang #if defined(PETSC_HAVE_ELEMENTAL) 1232cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatConvert_MPIDense_Elemental(Mat A, MatType newtype,MatReuse reuse,Mat *newmat) 12338baccfbdSHong Zhang { 12348ea901baSHong Zhang Mat mat_elemental; 12358ea901baSHong Zhang PetscErrorCode ierr; 123632d7a744SHong Zhang PetscScalar *v; 123732d7a744SHong Zhang PetscInt m=A->rmap->n,N=A->cmap->N,rstart=A->rmap->rstart,i,*rows,*cols; 12388ea901baSHong Zhang 12398baccfbdSHong Zhang PetscFunctionBegin; 1240378336b6SHong Zhang if (reuse == MAT_REUSE_MATRIX) { 1241378336b6SHong Zhang mat_elemental = *newmat; 1242378336b6SHong Zhang ierr = MatZeroEntries(*newmat);CHKERRQ(ierr); 1243378336b6SHong Zhang } else { 1244378336b6SHong Zhang ierr = MatCreate(PetscObjectComm((PetscObject)A), &mat_elemental);CHKERRQ(ierr); 1245378336b6SHong Zhang ierr = MatSetSizes(mat_elemental,PETSC_DECIDE,PETSC_DECIDE,A->rmap->N,A->cmap->N);CHKERRQ(ierr); 1246378336b6SHong Zhang ierr = MatSetType(mat_elemental,MATELEMENTAL);CHKERRQ(ierr); 1247378336b6SHong Zhang ierr = MatSetUp(mat_elemental);CHKERRQ(ierr); 124832d7a744SHong Zhang ierr = MatSetOption(mat_elemental,MAT_ROW_ORIENTED,PETSC_FALSE);CHKERRQ(ierr); 1249378336b6SHong Zhang } 1250378336b6SHong Zhang 125132d7a744SHong Zhang ierr = PetscMalloc2(m,&rows,N,&cols);CHKERRQ(ierr); 125232d7a744SHong Zhang for (i=0; i<N; i++) cols[i] = i; 125332d7a744SHong Zhang for (i=0; i<m; i++) rows[i] = rstart + i; 12548ea901baSHong Zhang 12558ea901baSHong Zhang /* PETSc-Elemental interaface uses axpy for setting off-processor entries, only ADD_VALUES is allowed */ 125632d7a744SHong Zhang ierr = MatDenseGetArray(A,&v);CHKERRQ(ierr); 125732d7a744SHong Zhang ierr = MatSetValues(mat_elemental,m,rows,N,cols,v,ADD_VALUES);CHKERRQ(ierr); 12588ea901baSHong Zhang ierr = MatAssemblyBegin(mat_elemental, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 12598ea901baSHong Zhang ierr = MatAssemblyEnd(mat_elemental, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 126032d7a744SHong Zhang ierr = MatDenseRestoreArray(A,&v);CHKERRQ(ierr); 126132d7a744SHong Zhang ierr = PetscFree2(rows,cols);CHKERRQ(ierr); 12628ea901baSHong Zhang 1263511c6705SHong Zhang if (reuse == MAT_INPLACE_MATRIX) { 126428be2f97SBarry Smith ierr = MatHeaderReplace(A,&mat_elemental);CHKERRQ(ierr); 12658ea901baSHong Zhang } else { 12668ea901baSHong Zhang *newmat = mat_elemental; 12678ea901baSHong Zhang } 12688baccfbdSHong Zhang PetscFunctionReturn(0); 12698baccfbdSHong Zhang } 127065b80a83SHong Zhang #endif 12718baccfbdSHong Zhang 1272af53bab2SHong Zhang static PetscErrorCode MatDenseGetColumn_MPIDense(Mat A,PetscInt col,PetscScalar **vals) 127386aefd0dSHong Zhang { 127486aefd0dSHong Zhang Mat_MPIDense *mat = (Mat_MPIDense*)A->data; 127586aefd0dSHong Zhang PetscErrorCode ierr; 127686aefd0dSHong Zhang 127786aefd0dSHong Zhang PetscFunctionBegin; 127886aefd0dSHong Zhang ierr = MatDenseGetColumn(mat->A,col,vals);CHKERRQ(ierr); 127986aefd0dSHong Zhang PetscFunctionReturn(0); 128086aefd0dSHong Zhang } 128186aefd0dSHong Zhang 1282af53bab2SHong Zhang static PetscErrorCode MatDenseRestoreColumn_MPIDense(Mat A,PetscScalar **vals) 128386aefd0dSHong Zhang { 128486aefd0dSHong Zhang Mat_MPIDense *mat = (Mat_MPIDense*)A->data; 128586aefd0dSHong Zhang PetscErrorCode ierr; 128686aefd0dSHong Zhang 128786aefd0dSHong Zhang PetscFunctionBegin; 128886aefd0dSHong Zhang ierr = MatDenseRestoreColumn(mat->A,vals);CHKERRQ(ierr); 128986aefd0dSHong Zhang PetscFunctionReturn(0); 129086aefd0dSHong Zhang } 129186aefd0dSHong Zhang 129294e2cb23SJakub Kruzik PetscErrorCode MatCreateMPIMatConcatenateSeqMat_MPIDense(MPI_Comm comm,Mat inmat,PetscInt n,MatReuse scall,Mat *outmat) 129394e2cb23SJakub Kruzik { 129494e2cb23SJakub Kruzik PetscErrorCode ierr; 129594e2cb23SJakub Kruzik Mat_MPIDense *mat; 129694e2cb23SJakub Kruzik PetscInt m,nloc,N; 129794e2cb23SJakub Kruzik 129894e2cb23SJakub Kruzik PetscFunctionBegin; 129994e2cb23SJakub Kruzik ierr = MatGetSize(inmat,&m,&N);CHKERRQ(ierr); 130094e2cb23SJakub Kruzik ierr = MatGetLocalSize(inmat,NULL,&nloc);CHKERRQ(ierr); 130194e2cb23SJakub Kruzik if (scall == MAT_INITIAL_MATRIX) { /* symbolic phase */ 130294e2cb23SJakub Kruzik PetscInt sum; 130394e2cb23SJakub Kruzik 130494e2cb23SJakub Kruzik if (n == PETSC_DECIDE) { 130594e2cb23SJakub Kruzik ierr = PetscSplitOwnership(comm,&n,&N);CHKERRQ(ierr); 130694e2cb23SJakub Kruzik } 130794e2cb23SJakub Kruzik /* Check sum(n) = N */ 130894e2cb23SJakub Kruzik ierr = MPIU_Allreduce(&n,&sum,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(ierr); 130994e2cb23SJakub Kruzik if (sum != N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Sum of local columns %D != global columns %D",sum,N); 131094e2cb23SJakub Kruzik 131194e2cb23SJakub Kruzik ierr = MatCreateDense(comm,m,n,PETSC_DETERMINE,N,NULL,outmat);CHKERRQ(ierr); 131294e2cb23SJakub Kruzik } 131394e2cb23SJakub Kruzik 131494e2cb23SJakub Kruzik /* numeric phase */ 131594e2cb23SJakub Kruzik mat = (Mat_MPIDense*)(*outmat)->data; 131694e2cb23SJakub Kruzik ierr = MatCopy(inmat,mat->A,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 131794e2cb23SJakub Kruzik ierr = MatAssemblyBegin(*outmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 131894e2cb23SJakub Kruzik ierr = MatAssemblyEnd(*outmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 131994e2cb23SJakub Kruzik PetscFunctionReturn(0); 132094e2cb23SJakub Kruzik } 132194e2cb23SJakub Kruzik 13228cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_MPIDense(Mat mat) 1323273d9f13SBarry Smith { 1324273d9f13SBarry Smith Mat_MPIDense *a; 1325dfbe8321SBarry Smith PetscErrorCode ierr; 1326273d9f13SBarry Smith 1327273d9f13SBarry Smith PetscFunctionBegin; 1328b00a9115SJed Brown ierr = PetscNewLog(mat,&a);CHKERRQ(ierr); 1329b0a32e0cSBarry Smith mat->data = (void*)a; 1330273d9f13SBarry Smith ierr = PetscMemcpy(mat->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 1331273d9f13SBarry Smith 1332273d9f13SBarry Smith mat->insertmode = NOT_SET_VALUES; 1333ce94432eSBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)mat),&a->rank);CHKERRQ(ierr); 1334ce94432eSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)mat),&a->size);CHKERRQ(ierr); 1335273d9f13SBarry Smith 1336273d9f13SBarry Smith /* build cache for off array entries formed */ 1337273d9f13SBarry Smith a->donotstash = PETSC_FALSE; 13382205254eSKarl Rupp 1339ce94432eSBarry Smith ierr = MatStashCreate_Private(PetscObjectComm((PetscObject)mat),1,&mat->stash);CHKERRQ(ierr); 1340273d9f13SBarry Smith 1341273d9f13SBarry Smith /* stuff used for matrix vector multiply */ 1342273d9f13SBarry Smith a->lvec = 0; 1343273d9f13SBarry Smith a->Mvctx = 0; 1344273d9f13SBarry Smith a->roworiented = PETSC_TRUE; 1345273d9f13SBarry Smith 134649a6ff4bSBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatDenseGetLDA_C",MatDenseGetLDA_MPIDense);CHKERRQ(ierr); 1347bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatDenseGetArray_C",MatDenseGetArray_MPIDense);CHKERRQ(ierr); 13488572280aSBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatDenseRestoreArray_C",MatDenseRestoreArray_MPIDense);CHKERRQ(ierr); 13498572280aSBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatDenseGetArrayRead_C",MatDenseGetArrayRead_MPIDense);CHKERRQ(ierr); 13508572280aSBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatDenseRestoreArrayRead_C",MatDenseRestoreArrayRead_MPIDense);CHKERRQ(ierr); 1351d3042a70SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatDensePlaceArray_C",MatDensePlaceArray_MPIDense);CHKERRQ(ierr); 1352d3042a70SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatDenseResetArray_C",MatDenseResetArray_MPIDense);CHKERRQ(ierr); 13538baccfbdSHong Zhang #if defined(PETSC_HAVE_ELEMENTAL) 13548baccfbdSHong Zhang ierr = PetscObjectComposeFunction((PetscObject)mat,"MatConvert_mpidense_elemental_C",MatConvert_MPIDense_Elemental);CHKERRQ(ierr); 13558baccfbdSHong Zhang #endif 1356bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMPIDenseSetPreallocation_C",MatMPIDenseSetPreallocation_MPIDense);CHKERRQ(ierr); 1357bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMatMult_mpiaij_mpidense_C",MatMatMult_MPIAIJ_MPIDense);CHKERRQ(ierr); 1358bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMatMultSymbolic_mpiaij_mpidense_C",MatMatMultSymbolic_MPIAIJ_MPIDense);CHKERRQ(ierr); 1359bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMatMultNumeric_mpiaij_mpidense_C",MatMatMultNumeric_MPIAIJ_MPIDense);CHKERRQ(ierr); 136052c5f739Sprj- ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMatMult_nest_mpidense_C",MatMatMult_Nest_Dense);CHKERRQ(ierr); 136152c5f739Sprj- ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMatMultSymbolic_nest_mpidense_C",MatMatMultSymbolic_Nest_Dense);CHKERRQ(ierr); 136252c5f739Sprj- ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMatMultNumeric_nest_mpidense_C",MatMatMultNumeric_Nest_Dense);CHKERRQ(ierr); 13638949adfdSHong Zhang 13648949adfdSHong Zhang ierr = PetscObjectComposeFunction((PetscObject)mat,"MatTransposeMatMult_mpiaij_mpidense_C",MatTransposeMatMult_MPIAIJ_MPIDense);CHKERRQ(ierr); 13658949adfdSHong Zhang ierr = PetscObjectComposeFunction((PetscObject)mat,"MatTransposeMatMultSymbolic_mpiaij_mpidense_C",MatTransposeMatMultSymbolic_MPIAIJ_MPIDense);CHKERRQ(ierr); 13668949adfdSHong Zhang ierr = PetscObjectComposeFunction((PetscObject)mat,"MatTransposeMatMultNumeric_mpiaij_mpidense_C",MatTransposeMatMultNumeric_MPIAIJ_MPIDense);CHKERRQ(ierr); 1367af53bab2SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)mat,"MatDenseGetColumn_C",MatDenseGetColumn_MPIDense);CHKERRQ(ierr); 1368af53bab2SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)mat,"MatDenseRestoreColumn_C",MatDenseRestoreColumn_MPIDense);CHKERRQ(ierr); 136938aed534SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)mat,MATMPIDENSE);CHKERRQ(ierr); 1370273d9f13SBarry Smith PetscFunctionReturn(0); 1371273d9f13SBarry Smith } 1372273d9f13SBarry Smith 1373209238afSKris Buschelman /*MC 1374002d173eSKris Buschelman MATDENSE - MATDENSE = "dense" - A matrix type to be used for dense matrices. 1375209238afSKris Buschelman 1376209238afSKris Buschelman This matrix type is identical to MATSEQDENSE when constructed with a single process communicator, 1377209238afSKris Buschelman and MATMPIDENSE otherwise. 1378209238afSKris Buschelman 1379209238afSKris Buschelman Options Database Keys: 1380209238afSKris Buschelman . -mat_type dense - sets the matrix type to "dense" during a call to MatSetFromOptions() 1381209238afSKris Buschelman 1382209238afSKris Buschelman Level: beginner 1383209238afSKris Buschelman 138401b82886SBarry Smith 138552c5f739Sprj- .seealso: MATSEQDENSE,MATMPIDENSE 1386209238afSKris Buschelman M*/ 1387209238afSKris Buschelman 1388273d9f13SBarry Smith /*@C 1389273d9f13SBarry Smith MatMPIDenseSetPreallocation - Sets the array used to store the matrix entries 1390273d9f13SBarry Smith 1391273d9f13SBarry Smith Not collective 1392273d9f13SBarry Smith 1393273d9f13SBarry Smith Input Parameters: 13941c4f3114SJed Brown . B - the matrix 13950298fd71SBarry Smith - data - optional location of matrix data. Set data=NULL for PETSc 1396273d9f13SBarry Smith to control all matrix memory allocation. 1397273d9f13SBarry Smith 1398273d9f13SBarry Smith Notes: 1399273d9f13SBarry Smith The dense format is fully compatible with standard Fortran 77 1400273d9f13SBarry Smith storage by columns. 1401273d9f13SBarry Smith 1402273d9f13SBarry Smith The data input variable is intended primarily for Fortran programmers 1403273d9f13SBarry Smith who wish to allocate their own matrix memory space. Most users should 14040298fd71SBarry Smith set data=NULL. 1405273d9f13SBarry Smith 1406273d9f13SBarry Smith Level: intermediate 1407273d9f13SBarry Smith 1408273d9f13SBarry Smith .seealso: MatCreate(), MatCreateSeqDense(), MatSetValues() 1409273d9f13SBarry Smith @*/ 14101c4f3114SJed Brown PetscErrorCode MatMPIDenseSetPreallocation(Mat B,PetscScalar *data) 1411273d9f13SBarry Smith { 14124ac538c5SBarry Smith PetscErrorCode ierr; 1413273d9f13SBarry Smith 1414273d9f13SBarry Smith PetscFunctionBegin; 14151c4f3114SJed Brown ierr = PetscTryMethod(B,"MatMPIDenseSetPreallocation_C",(Mat,PetscScalar*),(B,data));CHKERRQ(ierr); 1416273d9f13SBarry Smith PetscFunctionReturn(0); 1417273d9f13SBarry Smith } 1418273d9f13SBarry Smith 1419d3042a70SBarry Smith /*@ 1420d3042a70SBarry Smith MatDensePlaceArray - Allows one to replace the array in a dense array with an 1421d3042a70SBarry Smith array provided by the user. This is useful to avoid copying an array 1422d3042a70SBarry Smith into a matrix 1423d3042a70SBarry Smith 1424d3042a70SBarry Smith Not Collective 1425d3042a70SBarry Smith 1426d3042a70SBarry Smith Input Parameters: 1427d3042a70SBarry Smith + mat - the matrix 1428d3042a70SBarry Smith - array - the array in column major order 1429d3042a70SBarry Smith 1430d3042a70SBarry Smith Notes: 1431d3042a70SBarry 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 1432d3042a70SBarry Smith freed when the matrix is destroyed. 1433d3042a70SBarry Smith 1434d3042a70SBarry Smith Level: developer 1435d3042a70SBarry Smith 1436d3042a70SBarry Smith .seealso: MatDenseGetArray(), MatDenseResetArray(), VecPlaceArray(), VecGetArray(), VecRestoreArray(), VecReplaceArray(), VecResetArray() 1437d3042a70SBarry Smith 1438d3042a70SBarry Smith @*/ 1439d3042a70SBarry Smith PetscErrorCode MatDensePlaceArray(Mat mat,const PetscScalar array[]) 1440d3042a70SBarry Smith { 1441d3042a70SBarry Smith PetscErrorCode ierr; 1442d3042a70SBarry Smith PetscFunctionBegin; 1443d3042a70SBarry Smith ierr = PetscUseMethod(mat,"MatDensePlaceArray_C",(Mat,const PetscScalar*),(mat,array));CHKERRQ(ierr); 1444d3042a70SBarry Smith ierr = PetscObjectStateIncrease((PetscObject)mat);CHKERRQ(ierr); 1445d3042a70SBarry Smith PetscFunctionReturn(0); 1446d3042a70SBarry Smith } 1447d3042a70SBarry Smith 1448d3042a70SBarry Smith /*@ 1449d3042a70SBarry Smith MatDenseResetArray - Resets the matrix array to that it previously had before the call to MatDensePlaceArray() 1450d3042a70SBarry Smith 1451d3042a70SBarry Smith Not Collective 1452d3042a70SBarry Smith 1453d3042a70SBarry Smith Input Parameters: 1454d3042a70SBarry Smith . mat - the matrix 1455d3042a70SBarry Smith 1456d3042a70SBarry Smith Notes: 1457d3042a70SBarry Smith You can only call this after a call to MatDensePlaceArray() 1458d3042a70SBarry Smith 1459d3042a70SBarry Smith Level: developer 1460d3042a70SBarry Smith 1461d3042a70SBarry Smith .seealso: MatDenseGetArray(), MatDensePlaceArray(), VecPlaceArray(), VecGetArray(), VecRestoreArray(), VecReplaceArray(), VecResetArray() 1462d3042a70SBarry Smith 1463d3042a70SBarry Smith @*/ 1464d3042a70SBarry Smith PetscErrorCode MatDenseResetArray(Mat mat) 1465d3042a70SBarry Smith { 1466d3042a70SBarry Smith PetscErrorCode ierr; 1467d3042a70SBarry Smith PetscFunctionBegin; 1468d3042a70SBarry Smith ierr = PetscUseMethod(mat,"MatDenseResetArray_C",(Mat),(mat));CHKERRQ(ierr); 1469d3042a70SBarry Smith ierr = PetscObjectStateIncrease((PetscObject)mat);CHKERRQ(ierr); 1470d3042a70SBarry Smith PetscFunctionReturn(0); 1471d3042a70SBarry Smith } 1472d3042a70SBarry Smith 14738965ea79SLois Curfman McInnes /*@C 147469b1f4b7SBarry Smith MatCreateDense - Creates a parallel matrix in dense format. 14758965ea79SLois Curfman McInnes 1476d083f849SBarry Smith Collective 1477db81eaa0SLois Curfman McInnes 14788965ea79SLois Curfman McInnes Input Parameters: 1479db81eaa0SLois Curfman McInnes + comm - MPI communicator 14808965ea79SLois Curfman McInnes . m - number of local rows (or PETSC_DECIDE to have calculated if M is given) 1481db81eaa0SLois Curfman McInnes . n - number of local columns (or PETSC_DECIDE to have calculated if N is given) 14828965ea79SLois Curfman McInnes . M - number of global rows (or PETSC_DECIDE to have calculated if m is given) 1483db81eaa0SLois Curfman McInnes . N - number of global columns (or PETSC_DECIDE to have calculated if n is given) 14846cfe35ddSJose E. Roman - data - optional location of matrix data. Set data=NULL (PETSC_NULL_SCALAR for Fortran users) for PETSc 1485dfc5480cSLois Curfman McInnes to control all matrix memory allocation. 14868965ea79SLois Curfman McInnes 14878965ea79SLois Curfman McInnes Output Parameter: 1488477f1c0bSLois Curfman McInnes . A - the matrix 14898965ea79SLois Curfman McInnes 1490b259b22eSLois Curfman McInnes Notes: 149139ddd567SLois Curfman McInnes The dense format is fully compatible with standard Fortran 77 149239ddd567SLois Curfman McInnes storage by columns. 14938965ea79SLois Curfman McInnes 149418f449edSLois Curfman McInnes The data input variable is intended primarily for Fortran programmers 149518f449edSLois Curfman McInnes who wish to allocate their own matrix memory space. Most users should 14966cfe35ddSJose E. Roman set data=NULL (PETSC_NULL_SCALAR for Fortran users). 149718f449edSLois Curfman McInnes 14988965ea79SLois Curfman McInnes The user MUST specify either the local or global matrix dimensions 14998965ea79SLois Curfman McInnes (possibly both). 15008965ea79SLois Curfman McInnes 1501027ccd11SLois Curfman McInnes Level: intermediate 1502027ccd11SLois Curfman McInnes 150339ddd567SLois Curfman McInnes .seealso: MatCreate(), MatCreateSeqDense(), MatSetValues() 15048965ea79SLois Curfman McInnes @*/ 150569b1f4b7SBarry Smith PetscErrorCode MatCreateDense(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscScalar *data,Mat *A) 15068965ea79SLois Curfman McInnes { 15076849ba73SBarry Smith PetscErrorCode ierr; 150813f74950SBarry Smith PetscMPIInt size; 15098965ea79SLois Curfman McInnes 15103a40ed3dSBarry Smith PetscFunctionBegin; 1511f69a0ea3SMatthew Knepley ierr = MatCreate(comm,A);CHKERRQ(ierr); 1512*8491ab44SLisandro Dalcin PetscValidLogicalCollectiveBool(*A,!!data,6); 1513f69a0ea3SMatthew Knepley ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr); 1514273d9f13SBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1515273d9f13SBarry Smith if (size > 1) { 1516273d9f13SBarry Smith ierr = MatSetType(*A,MATMPIDENSE);CHKERRQ(ierr); 1517273d9f13SBarry Smith ierr = MatMPIDenseSetPreallocation(*A,data);CHKERRQ(ierr); 15186cfe35ddSJose E. Roman if (data) { /* user provided data array, so no need to assemble */ 15196cfe35ddSJose E. Roman ierr = MatSetUpMultiply_MPIDense(*A);CHKERRQ(ierr); 15206cfe35ddSJose E. Roman (*A)->assembled = PETSC_TRUE; 15216cfe35ddSJose E. Roman } 1522273d9f13SBarry Smith } else { 1523273d9f13SBarry Smith ierr = MatSetType(*A,MATSEQDENSE);CHKERRQ(ierr); 1524273d9f13SBarry Smith ierr = MatSeqDenseSetPreallocation(*A,data);CHKERRQ(ierr); 15258c469469SLois Curfman McInnes } 15263a40ed3dSBarry Smith PetscFunctionReturn(0); 15278965ea79SLois Curfman McInnes } 15288965ea79SLois Curfman McInnes 15296849ba73SBarry Smith static PetscErrorCode MatDuplicate_MPIDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat) 15308965ea79SLois Curfman McInnes { 15318965ea79SLois Curfman McInnes Mat mat; 15323501a2bdSLois Curfman McInnes Mat_MPIDense *a,*oldmat = (Mat_MPIDense*)A->data; 1533dfbe8321SBarry Smith PetscErrorCode ierr; 15348965ea79SLois Curfman McInnes 15353a40ed3dSBarry Smith PetscFunctionBegin; 15368965ea79SLois Curfman McInnes *newmat = 0; 1537ce94432eSBarry Smith ierr = MatCreate(PetscObjectComm((PetscObject)A),&mat);CHKERRQ(ierr); 1538d0f46423SBarry Smith ierr = MatSetSizes(mat,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr); 15397adad957SLisandro Dalcin ierr = MatSetType(mat,((PetscObject)A)->type_name);CHKERRQ(ierr); 1540834f8fabSBarry Smith a = (Mat_MPIDense*)mat->data; 15415aa7edbeSHong Zhang 1542d5f3da31SBarry Smith mat->factortype = A->factortype; 1543c456f294SBarry Smith mat->assembled = PETSC_TRUE; 1544273d9f13SBarry Smith mat->preallocated = PETSC_TRUE; 15458965ea79SLois Curfman McInnes 15468965ea79SLois Curfman McInnes a->size = oldmat->size; 15478965ea79SLois Curfman McInnes a->rank = oldmat->rank; 1548e0fa3b82SLois Curfman McInnes mat->insertmode = NOT_SET_VALUES; 1549b9b97703SBarry Smith a->nvec = oldmat->nvec; 15503782ba37SSatish Balay a->donotstash = oldmat->donotstash; 1551e04c1aa4SHong Zhang 15521e1e43feSBarry Smith ierr = PetscLayoutReference(A->rmap,&mat->rmap);CHKERRQ(ierr); 15531e1e43feSBarry Smith ierr = PetscLayoutReference(A->cmap,&mat->cmap);CHKERRQ(ierr); 15548965ea79SLois Curfman McInnes 1555329f5518SBarry Smith ierr = MatSetUpMultiply_MPIDense(mat);CHKERRQ(ierr); 15565609ef8eSBarry Smith ierr = MatDuplicate(oldmat->A,cpvalues,&a->A);CHKERRQ(ierr); 15573bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)a->A);CHKERRQ(ierr); 155801b82886SBarry Smith 15598965ea79SLois Curfman McInnes *newmat = mat; 15603a40ed3dSBarry Smith PetscFunctionReturn(0); 15618965ea79SLois Curfman McInnes } 15628965ea79SLois Curfman McInnes 1563eb91f321SVaclav Hapla PetscErrorCode MatLoad_MPIDense(Mat newMat, PetscViewer viewer) 1564eb91f321SVaclav Hapla { 1565eb91f321SVaclav Hapla PetscErrorCode ierr; 156687d5ce66SSatish Balay PetscBool isbinary; 156787d5ce66SSatish Balay #if defined(PETSC_HAVE_HDF5) 156887d5ce66SSatish Balay PetscBool ishdf5; 156987d5ce66SSatish Balay #endif 1570eb91f321SVaclav Hapla 1571eb91f321SVaclav Hapla PetscFunctionBegin; 1572eb91f321SVaclav Hapla PetscValidHeaderSpecific(newMat,MAT_CLASSID,1); 1573eb91f321SVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 1574eb91f321SVaclav Hapla /* force binary viewer to load .info file if it has not yet done so */ 1575eb91f321SVaclav Hapla ierr = PetscViewerSetUp(viewer);CHKERRQ(ierr); 1576eb91f321SVaclav Hapla ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr); 157787d5ce66SSatish Balay #if defined(PETSC_HAVE_HDF5) 1578eb91f321SVaclav Hapla ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERHDF5, &ishdf5);CHKERRQ(ierr); 157987d5ce66SSatish Balay #endif 1580eb91f321SVaclav Hapla if (isbinary) { 1581*8491ab44SLisandro Dalcin ierr = MatLoad_Dense_Binary(newMat,viewer);CHKERRQ(ierr); 1582eb91f321SVaclav Hapla #if defined(PETSC_HAVE_HDF5) 158387d5ce66SSatish Balay } else if (ishdf5) { 1584eb91f321SVaclav Hapla ierr = MatLoad_Dense_HDF5(newMat,viewer);CHKERRQ(ierr); 1585eb91f321SVaclav Hapla #endif 158687d5ce66SSatish Balay } else SETERRQ2(PetscObjectComm((PetscObject)newMat),PETSC_ERR_SUP,"Viewer type %s not yet supported for reading %s matrices",((PetscObject)viewer)->type_name,((PetscObject)newMat)->type_name); 1587eb91f321SVaclav Hapla PetscFunctionReturn(0); 1588eb91f321SVaclav Hapla } 1589eb91f321SVaclav Hapla 1590ace3abfcSBarry Smith PetscErrorCode MatEqual_MPIDense(Mat A,Mat B,PetscBool *flag) 15916e4ee0c6SHong Zhang { 15926e4ee0c6SHong Zhang Mat_MPIDense *matB = (Mat_MPIDense*)B->data,*matA = (Mat_MPIDense*)A->data; 15936e4ee0c6SHong Zhang Mat a,b; 1594ace3abfcSBarry Smith PetscBool flg; 15956e4ee0c6SHong Zhang PetscErrorCode ierr; 159690ace30eSBarry Smith 15976e4ee0c6SHong Zhang PetscFunctionBegin; 15986e4ee0c6SHong Zhang a = matA->A; 15996e4ee0c6SHong Zhang b = matB->A; 16006e4ee0c6SHong Zhang ierr = MatEqual(a,b,&flg);CHKERRQ(ierr); 1601b2566f29SBarry Smith ierr = MPIU_Allreduce(&flg,flag,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 16026e4ee0c6SHong Zhang PetscFunctionReturn(0); 16036e4ee0c6SHong Zhang } 160490ace30eSBarry Smith 1605baa3c1c6SHong Zhang PetscErrorCode MatDestroy_MatTransMatMult_MPIDense_MPIDense(Mat A) 1606baa3c1c6SHong Zhang { 1607baa3c1c6SHong Zhang PetscErrorCode ierr; 1608baa3c1c6SHong Zhang Mat_MPIDense *a = (Mat_MPIDense*)A->data; 1609baa3c1c6SHong Zhang Mat_TransMatMultDense *atb = a->atbdense; 1610baa3c1c6SHong Zhang 1611baa3c1c6SHong Zhang PetscFunctionBegin; 1612c5ef1628SHong Zhang ierr = PetscFree3(atb->sendbuf,atb->atbarray,atb->recvcounts);CHKERRQ(ierr); 1613baa3c1c6SHong Zhang ierr = (atb->destroy)(A);CHKERRQ(ierr); 1614baa3c1c6SHong Zhang ierr = PetscFree(atb);CHKERRQ(ierr); 1615baa3c1c6SHong Zhang PetscFunctionReturn(0); 1616baa3c1c6SHong Zhang } 1617baa3c1c6SHong Zhang 1618cc48ffa7SToby Isaac PetscErrorCode MatDestroy_MatMatTransMult_MPIDense_MPIDense(Mat A) 1619cc48ffa7SToby Isaac { 1620cc48ffa7SToby Isaac PetscErrorCode ierr; 1621cc48ffa7SToby Isaac Mat_MPIDense *a = (Mat_MPIDense*)A->data; 1622cc48ffa7SToby Isaac Mat_MatTransMultDense *abt = a->abtdense; 1623cc48ffa7SToby Isaac 1624cc48ffa7SToby Isaac PetscFunctionBegin; 1625cc48ffa7SToby Isaac ierr = PetscFree2(abt->buf[0],abt->buf[1]);CHKERRQ(ierr); 1626faa55883SToby Isaac ierr = PetscFree2(abt->recvcounts,abt->recvdispls);CHKERRQ(ierr); 1627cc48ffa7SToby Isaac ierr = (abt->destroy)(A);CHKERRQ(ierr); 1628cc48ffa7SToby Isaac ierr = PetscFree(abt);CHKERRQ(ierr); 1629cc48ffa7SToby Isaac PetscFunctionReturn(0); 1630cc48ffa7SToby Isaac } 1631cc48ffa7SToby Isaac 1632cb20be35SHong Zhang PetscErrorCode MatTransposeMatMultNumeric_MPIDense_MPIDense(Mat A,Mat B,Mat C) 1633cb20be35SHong Zhang { 1634baa3c1c6SHong Zhang Mat_MPIDense *a=(Mat_MPIDense*)A->data, *b=(Mat_MPIDense*)B->data, *c=(Mat_MPIDense*)C->data; 1635660d5466SHong Zhang Mat_SeqDense *aseq=(Mat_SeqDense*)(a->A)->data, *bseq=(Mat_SeqDense*)(b->A)->data; 1636baa3c1c6SHong Zhang Mat_TransMatMultDense *atb = c->atbdense; 1637cb20be35SHong Zhang PetscErrorCode ierr; 1638cb20be35SHong Zhang MPI_Comm comm; 1639d5017740SHong Zhang PetscMPIInt rank,size,*recvcounts=atb->recvcounts; 1640c5ef1628SHong Zhang PetscScalar *carray,*atbarray=atb->atbarray,*sendbuf=atb->sendbuf; 1641d5017740SHong Zhang PetscInt i,cN=C->cmap->N,cM=C->rmap->N,proc,k,j; 1642660d5466SHong Zhang PetscScalar _DOne=1.0,_DZero=0.0; 1643adc7a786SHong Zhang PetscBLASInt am,an,bn,aN; 1644e68c0b26SHong Zhang const PetscInt *ranges; 1645cb20be35SHong Zhang 1646cb20be35SHong Zhang PetscFunctionBegin; 1647cb20be35SHong Zhang ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 1648cb20be35SHong Zhang ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 1649cb20be35SHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1650e68c0b26SHong Zhang 1651c5ef1628SHong Zhang /* compute atbarray = aseq^T * bseq */ 1652660d5466SHong Zhang ierr = PetscBLASIntCast(a->A->cmap->n,&an);CHKERRQ(ierr); 1653660d5466SHong Zhang ierr = PetscBLASIntCast(b->A->cmap->n,&bn);CHKERRQ(ierr); 1654660d5466SHong Zhang ierr = PetscBLASIntCast(a->A->rmap->n,&am);CHKERRQ(ierr); 1655adc7a786SHong Zhang ierr = PetscBLASIntCast(A->cmap->N,&aN);CHKERRQ(ierr); 1656adc7a786SHong Zhang PetscStackCallBLAS("BLASgemm",BLASgemm_("T","N",&an,&bn,&am,&_DOne,aseq->v,&aseq->lda,bseq->v,&bseq->lda,&_DZero,atbarray,&aN)); 1657cb20be35SHong Zhang 1658cb20be35SHong Zhang ierr = MatGetOwnershipRanges(C,&ranges);CHKERRQ(ierr); 1659c5ef1628SHong Zhang for (i=0; i<size; i++) recvcounts[i] = (ranges[i+1] - ranges[i])*cN; 1660cb20be35SHong Zhang 1661660d5466SHong Zhang /* arrange atbarray into sendbuf */ 1662cb20be35SHong Zhang k = 0; 1663cb20be35SHong Zhang for (proc=0; proc<size; proc++) { 1664baa3c1c6SHong Zhang for (j=0; j<cN; j++) { 1665c5ef1628SHong Zhang for (i=ranges[proc]; i<ranges[proc+1]; i++) sendbuf[k++] = atbarray[i+j*cM]; 1666cb20be35SHong Zhang } 1667cb20be35SHong Zhang } 1668c5ef1628SHong Zhang /* sum all atbarray to local values of C */ 1669660d5466SHong Zhang ierr = MatDenseGetArray(c->A,&carray);CHKERRQ(ierr); 16703462b7efSHong Zhang ierr = MPI_Reduce_scatter(sendbuf,carray,recvcounts,MPIU_SCALAR,MPIU_SUM,comm);CHKERRQ(ierr); 1671660d5466SHong Zhang ierr = MatDenseRestoreArray(c->A,&carray);CHKERRQ(ierr); 1672cb20be35SHong Zhang PetscFunctionReturn(0); 1673cb20be35SHong Zhang } 1674cb20be35SHong Zhang 1675cb20be35SHong Zhang PetscErrorCode MatTransposeMatMultSymbolic_MPIDense_MPIDense(Mat A,Mat B,PetscReal fill,Mat *C) 1676cb20be35SHong Zhang { 1677cb20be35SHong Zhang PetscErrorCode ierr; 1678baa3c1c6SHong Zhang Mat Cdense; 1679cb20be35SHong Zhang MPI_Comm comm; 1680baa3c1c6SHong Zhang PetscMPIInt size; 1681660d5466SHong Zhang PetscInt cm=A->cmap->n,cM,cN=B->cmap->N; 1682baa3c1c6SHong Zhang Mat_MPIDense *c; 1683baa3c1c6SHong Zhang Mat_TransMatMultDense *atb; 1684cb20be35SHong Zhang 1685cb20be35SHong Zhang PetscFunctionBegin; 1686baa3c1c6SHong Zhang ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 1687cb20be35SHong Zhang if (A->rmap->rstart != B->rmap->rstart || A->rmap->rend != B->rmap->rend) { 1688cb20be35SHong 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); 1689cb20be35SHong Zhang } 1690cb20be35SHong Zhang 1691baa3c1c6SHong Zhang /* create matrix product Cdense */ 1692baa3c1c6SHong Zhang ierr = MatCreate(comm,&Cdense);CHKERRQ(ierr); 1693660d5466SHong Zhang ierr = MatSetSizes(Cdense,cm,B->cmap->n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr); 1694baa3c1c6SHong Zhang ierr = MatSetType(Cdense,MATMPIDENSE);CHKERRQ(ierr); 1695baa3c1c6SHong Zhang ierr = MatMPIDenseSetPreallocation(Cdense,NULL);CHKERRQ(ierr); 1696baa3c1c6SHong Zhang ierr = MatAssemblyBegin(Cdense,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1697baa3c1c6SHong Zhang ierr = MatAssemblyEnd(Cdense,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1698baa3c1c6SHong Zhang *C = Cdense; 1699baa3c1c6SHong Zhang 1700baa3c1c6SHong Zhang /* create data structure for reuse Cdense */ 1701baa3c1c6SHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1702baa3c1c6SHong Zhang ierr = PetscNew(&atb);CHKERRQ(ierr); 1703660d5466SHong Zhang cM = Cdense->rmap->N; 1704c5ef1628SHong Zhang ierr = PetscMalloc3(cM*cN,&atb->sendbuf,cM*cN,&atb->atbarray,size,&atb->recvcounts);CHKERRQ(ierr); 1705baa3c1c6SHong Zhang 1706baa3c1c6SHong Zhang c = (Mat_MPIDense*)Cdense->data; 1707baa3c1c6SHong Zhang c->atbdense = atb; 1708baa3c1c6SHong Zhang atb->destroy = Cdense->ops->destroy; 1709baa3c1c6SHong Zhang Cdense->ops->destroy = MatDestroy_MatTransMatMult_MPIDense_MPIDense; 1710cb20be35SHong Zhang PetscFunctionReturn(0); 1711cb20be35SHong Zhang } 1712cb20be35SHong Zhang 1713cb20be35SHong Zhang PetscErrorCode MatTransposeMatMult_MPIDense_MPIDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 1714cb20be35SHong Zhang { 1715cb20be35SHong Zhang PetscErrorCode ierr; 1716cb20be35SHong Zhang 1717cb20be35SHong Zhang PetscFunctionBegin; 1718cb20be35SHong Zhang if (scall == MAT_INITIAL_MATRIX) { 1719cb20be35SHong Zhang ierr = MatTransposeMatMultSymbolic_MPIDense_MPIDense(A,B,fill,C);CHKERRQ(ierr); 1720cb20be35SHong Zhang } 1721cb20be35SHong Zhang ierr = MatTransposeMatMultNumeric_MPIDense_MPIDense(A,B,*C);CHKERRQ(ierr); 1722cb20be35SHong Zhang PetscFunctionReturn(0); 1723cb20be35SHong Zhang } 1724320f2790SHong Zhang 1725cc48ffa7SToby Isaac static PetscErrorCode MatMatTransposeMultSymbolic_MPIDense_MPIDense(Mat A, Mat B, PetscReal fill, Mat *C) 1726cc48ffa7SToby Isaac { 1727cc48ffa7SToby Isaac PetscErrorCode ierr; 1728cc48ffa7SToby Isaac Mat Cdense; 1729cc48ffa7SToby Isaac MPI_Comm comm; 1730cc48ffa7SToby Isaac PetscMPIInt i, size; 1731cc48ffa7SToby Isaac PetscInt maxRows, bufsiz; 1732cc48ffa7SToby Isaac Mat_MPIDense *c; 1733cc48ffa7SToby Isaac PetscMPIInt tag; 1734faa55883SToby Isaac const char *algTypes[2] = {"allgatherv","cyclic"}; 1735faa55883SToby Isaac PetscInt alg, nalg = 2; 1736cc48ffa7SToby Isaac Mat_MatTransMultDense *abt; 1737cc48ffa7SToby Isaac 1738cc48ffa7SToby Isaac PetscFunctionBegin; 1739cc48ffa7SToby Isaac ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 1740faa55883SToby Isaac alg = 0; /* default is allgatherv */ 1741faa55883SToby Isaac ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"MatMatTransposeMult","Mat");CHKERRQ(ierr); 1742faa55883SToby Isaac ierr = PetscOptionsEList("-matmattransmult_mpidense_mpidense_via","Algorithmic approach","MatMatTransposeMult",algTypes,nalg,algTypes[0],&alg,NULL);CHKERRQ(ierr); 1743faa55883SToby Isaac ierr = PetscOptionsEnd();CHKERRQ(ierr); 1744cc48ffa7SToby Isaac if (A->cmap->N != B->cmap->N) { 1745cc48ffa7SToby Isaac SETERRQ2(comm,PETSC_ERR_ARG_SIZ,"Matrix global column dimensions are incompatible, A (%D) != B (%D)",A->cmap->N,B->cmap->N); 1746cc48ffa7SToby Isaac } 1747cc48ffa7SToby Isaac 1748cc48ffa7SToby Isaac /* create matrix product Cdense */ 1749cc48ffa7SToby Isaac ierr = MatCreate(comm,&Cdense);CHKERRQ(ierr); 1750cc48ffa7SToby Isaac ierr = MatSetSizes(Cdense,A->rmap->n,B->rmap->n,A->rmap->N,B->rmap->N);CHKERRQ(ierr); 1751cc48ffa7SToby Isaac ierr = MatSetType(Cdense,MATMPIDENSE);CHKERRQ(ierr); 1752cc48ffa7SToby Isaac ierr = MatMPIDenseSetPreallocation(Cdense,NULL);CHKERRQ(ierr); 1753cc48ffa7SToby Isaac ierr = PetscObjectGetNewTag((PetscObject)Cdense, &tag);CHKERRQ(ierr); 1754cc48ffa7SToby Isaac ierr = MatAssemblyBegin(Cdense,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1755cc48ffa7SToby Isaac ierr = MatAssemblyEnd(Cdense,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1756cc48ffa7SToby Isaac *C = Cdense; 1757cc48ffa7SToby Isaac 1758cc48ffa7SToby Isaac /* create data structure for reuse Cdense */ 1759cc48ffa7SToby Isaac ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1760cc48ffa7SToby Isaac ierr = PetscNew(&abt);CHKERRQ(ierr); 1761cc48ffa7SToby Isaac abt->tag = tag; 1762faa55883SToby Isaac abt->alg = alg; 1763faa55883SToby Isaac switch (alg) { 1764faa55883SToby Isaac case 1: 1765cc48ffa7SToby Isaac for (maxRows = 0, i = 0; i < size; i++) maxRows = PetscMax(maxRows, (B->rmap->range[i + 1] - B->rmap->range[i])); 1766cc48ffa7SToby Isaac bufsiz = A->cmap->N * maxRows; 1767cc48ffa7SToby Isaac ierr = PetscMalloc2(bufsiz,&(abt->buf[0]),bufsiz,&(abt->buf[1]));CHKERRQ(ierr); 1768faa55883SToby Isaac break; 1769faa55883SToby Isaac default: 1770faa55883SToby Isaac ierr = PetscMalloc2(B->rmap->n * B->cmap->N, &(abt->buf[0]), B->rmap->N * B->cmap->N, &(abt->buf[1]));CHKERRQ(ierr); 1771faa55883SToby Isaac ierr = PetscMalloc2(size,&(abt->recvcounts),size+1,&(abt->recvdispls));CHKERRQ(ierr); 1772faa55883SToby Isaac for (i = 0; i <= size; i++) abt->recvdispls[i] = B->rmap->range[i] * A->cmap->N; 1773faa55883SToby Isaac for (i = 0; i < size; i++) abt->recvcounts[i] = abt->recvdispls[i + 1] - abt->recvdispls[i]; 1774faa55883SToby Isaac break; 1775faa55883SToby Isaac } 1776cc48ffa7SToby Isaac 1777cc48ffa7SToby Isaac c = (Mat_MPIDense*)Cdense->data; 1778cc48ffa7SToby Isaac c->abtdense = abt; 1779cc48ffa7SToby Isaac abt->destroy = Cdense->ops->destroy; 1780cc48ffa7SToby Isaac Cdense->ops->destroy = MatDestroy_MatMatTransMult_MPIDense_MPIDense; 1781cc48ffa7SToby Isaac PetscFunctionReturn(0); 1782cc48ffa7SToby Isaac } 1783cc48ffa7SToby Isaac 1784faa55883SToby Isaac static PetscErrorCode MatMatTransposeMultNumeric_MPIDense_MPIDense_Cyclic(Mat A, Mat B, Mat C) 1785cc48ffa7SToby Isaac { 1786cc48ffa7SToby Isaac Mat_MPIDense *a=(Mat_MPIDense*)A->data, *b=(Mat_MPIDense*)B->data, *c=(Mat_MPIDense*)C->data; 1787cc48ffa7SToby Isaac Mat_SeqDense *aseq=(Mat_SeqDense*)(a->A)->data, *bseq=(Mat_SeqDense*)(b->A)->data; 1788cc48ffa7SToby Isaac Mat_SeqDense *cseq=(Mat_SeqDense*)(c->A)->data; 1789cc48ffa7SToby Isaac Mat_MatTransMultDense *abt = c->abtdense; 1790cc48ffa7SToby Isaac PetscErrorCode ierr; 1791cc48ffa7SToby Isaac MPI_Comm comm; 1792cc48ffa7SToby Isaac PetscMPIInt rank,size, sendsiz, recvsiz, sendto, recvfrom, recvisfrom; 1793ce36caf3SSatish Balay PetscScalar *sendbuf, *recvbuf=0, *carray; 1794cc48ffa7SToby Isaac PetscInt i,cK=A->cmap->N,k,j,bn; 1795cc48ffa7SToby Isaac PetscScalar _DOne=1.0,_DZero=0.0; 1796cc48ffa7SToby Isaac PetscBLASInt cm, cn, ck; 1797cc48ffa7SToby Isaac MPI_Request reqs[2]; 1798cc48ffa7SToby Isaac const PetscInt *ranges; 1799cc48ffa7SToby Isaac 1800cc48ffa7SToby Isaac PetscFunctionBegin; 1801cc48ffa7SToby Isaac ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 1802cc48ffa7SToby Isaac ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 1803cc48ffa7SToby Isaac ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1804cc48ffa7SToby Isaac 1805cc48ffa7SToby Isaac ierr = MatGetOwnershipRanges(B,&ranges);CHKERRQ(ierr); 1806cc48ffa7SToby Isaac bn = B->rmap->n; 1807cc48ffa7SToby Isaac if (bseq->lda == bn) { 1808cc48ffa7SToby Isaac sendbuf = bseq->v; 1809cc48ffa7SToby Isaac } else { 1810cc48ffa7SToby Isaac sendbuf = abt->buf[0]; 1811cc48ffa7SToby Isaac for (k = 0, i = 0; i < cK; i++) { 1812cc48ffa7SToby Isaac for (j = 0; j < bn; j++, k++) { 1813cc48ffa7SToby Isaac sendbuf[k] = bseq->v[i * bseq->lda + j]; 1814cc48ffa7SToby Isaac } 1815cc48ffa7SToby Isaac } 1816cc48ffa7SToby Isaac } 1817cc48ffa7SToby Isaac if (size > 1) { 1818cc48ffa7SToby Isaac sendto = (rank + size - 1) % size; 1819cc48ffa7SToby Isaac recvfrom = (rank + size + 1) % size; 1820cc48ffa7SToby Isaac } else { 1821cc48ffa7SToby Isaac sendto = recvfrom = 0; 1822cc48ffa7SToby Isaac } 1823cc48ffa7SToby Isaac ierr = PetscBLASIntCast(cK,&ck);CHKERRQ(ierr); 1824cc48ffa7SToby Isaac ierr = PetscBLASIntCast(c->A->rmap->n,&cm);CHKERRQ(ierr); 1825cc48ffa7SToby Isaac recvisfrom = rank; 1826cc48ffa7SToby Isaac for (i = 0; i < size; i++) { 1827cc48ffa7SToby Isaac /* we have finished receiving in sending, bufs can be read/modified */ 1828cc48ffa7SToby Isaac PetscInt nextrecvisfrom = (recvisfrom + 1) % size; /* which process the next recvbuf will originate on */ 1829cc48ffa7SToby Isaac PetscInt nextbn = ranges[nextrecvisfrom + 1] - ranges[nextrecvisfrom]; 1830cc48ffa7SToby Isaac 1831cc48ffa7SToby Isaac if (nextrecvisfrom != rank) { 1832cc48ffa7SToby Isaac /* start the cyclic sends from sendbuf, to recvbuf (which will switch to sendbuf) */ 1833cc48ffa7SToby Isaac sendsiz = cK * bn; 1834cc48ffa7SToby Isaac recvsiz = cK * nextbn; 1835cc48ffa7SToby Isaac recvbuf = (i & 1) ? abt->buf[0] : abt->buf[1]; 1836cc48ffa7SToby Isaac ierr = MPI_Isend(sendbuf, sendsiz, MPIU_SCALAR, sendto, abt->tag, comm, &reqs[0]);CHKERRQ(ierr); 1837cc48ffa7SToby Isaac ierr = MPI_Irecv(recvbuf, recvsiz, MPIU_SCALAR, recvfrom, abt->tag, comm, &reqs[1]);CHKERRQ(ierr); 1838cc48ffa7SToby Isaac } 1839cc48ffa7SToby Isaac 1840cc48ffa7SToby Isaac /* local aseq * sendbuf^T */ 1841cc48ffa7SToby Isaac ierr = PetscBLASIntCast(ranges[recvisfrom + 1] - ranges[recvisfrom], &cn);CHKERRQ(ierr); 1842cc48ffa7SToby Isaac carray = &cseq->v[cseq->lda * ranges[recvisfrom]]; 1843cc48ffa7SToby Isaac PetscStackCallBLAS("BLASgemm",BLASgemm_("N","T",&cm,&cn,&ck,&_DOne,aseq->v,&aseq->lda,sendbuf,&cn,&_DZero,carray,&cseq->lda)); 1844cc48ffa7SToby Isaac 1845cc48ffa7SToby Isaac if (nextrecvisfrom != rank) { 1846cc48ffa7SToby Isaac /* wait for the sends and receives to complete, swap sendbuf and recvbuf */ 1847cc48ffa7SToby Isaac ierr = MPI_Waitall(2, reqs, MPI_STATUSES_IGNORE);CHKERRQ(ierr); 1848cc48ffa7SToby Isaac } 1849cc48ffa7SToby Isaac bn = nextbn; 1850cc48ffa7SToby Isaac recvisfrom = nextrecvisfrom; 1851cc48ffa7SToby Isaac sendbuf = recvbuf; 1852cc48ffa7SToby Isaac } 1853cc48ffa7SToby Isaac PetscFunctionReturn(0); 1854cc48ffa7SToby Isaac } 1855cc48ffa7SToby Isaac 1856faa55883SToby Isaac static PetscErrorCode MatMatTransposeMultNumeric_MPIDense_MPIDense_Allgatherv(Mat A, Mat B, Mat C) 1857faa55883SToby Isaac { 1858faa55883SToby Isaac Mat_MPIDense *a=(Mat_MPIDense*)A->data, *b=(Mat_MPIDense*)B->data, *c=(Mat_MPIDense*)C->data; 1859faa55883SToby Isaac Mat_SeqDense *aseq=(Mat_SeqDense*)(a->A)->data, *bseq=(Mat_SeqDense*)(b->A)->data; 1860faa55883SToby Isaac Mat_SeqDense *cseq=(Mat_SeqDense*)(c->A)->data; 1861faa55883SToby Isaac Mat_MatTransMultDense *abt = c->abtdense; 1862faa55883SToby Isaac PetscErrorCode ierr; 1863faa55883SToby Isaac MPI_Comm comm; 1864faa55883SToby Isaac PetscMPIInt rank,size; 1865faa55883SToby Isaac PetscScalar *sendbuf, *recvbuf; 1866faa55883SToby Isaac PetscInt i,cK=A->cmap->N,k,j,bn; 1867faa55883SToby Isaac PetscScalar _DOne=1.0,_DZero=0.0; 1868faa55883SToby Isaac PetscBLASInt cm, cn, ck; 1869faa55883SToby Isaac 1870faa55883SToby Isaac PetscFunctionBegin; 1871faa55883SToby Isaac ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 1872faa55883SToby Isaac ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 1873faa55883SToby Isaac ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1874faa55883SToby Isaac 1875faa55883SToby Isaac /* copy transpose of B into buf[0] */ 1876faa55883SToby Isaac bn = B->rmap->n; 1877faa55883SToby Isaac sendbuf = abt->buf[0]; 1878faa55883SToby Isaac recvbuf = abt->buf[1]; 1879faa55883SToby Isaac for (k = 0, j = 0; j < bn; j++) { 1880faa55883SToby Isaac for (i = 0; i < cK; i++, k++) { 1881faa55883SToby Isaac sendbuf[k] = bseq->v[i * bseq->lda + j]; 1882faa55883SToby Isaac } 1883faa55883SToby Isaac } 1884faa55883SToby Isaac ierr = MPI_Allgatherv(sendbuf, bn * cK, MPIU_SCALAR, recvbuf, abt->recvcounts, abt->recvdispls, MPIU_SCALAR, comm);CHKERRQ(ierr); 1885faa55883SToby Isaac ierr = PetscBLASIntCast(cK,&ck);CHKERRQ(ierr); 1886faa55883SToby Isaac ierr = PetscBLASIntCast(c->A->rmap->n,&cm);CHKERRQ(ierr); 1887faa55883SToby Isaac ierr = PetscBLASIntCast(c->A->cmap->n,&cn);CHKERRQ(ierr); 1888faa55883SToby Isaac PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&cm,&cn,&ck,&_DOne,aseq->v,&aseq->lda,recvbuf,&ck,&_DZero,cseq->v,&cseq->lda)); 1889faa55883SToby Isaac PetscFunctionReturn(0); 1890faa55883SToby Isaac } 1891faa55883SToby Isaac 1892faa55883SToby Isaac static PetscErrorCode MatMatTransposeMultNumeric_MPIDense_MPIDense(Mat A, Mat B, Mat C) 1893faa55883SToby Isaac { 1894faa55883SToby Isaac Mat_MPIDense *c=(Mat_MPIDense*)C->data; 1895faa55883SToby Isaac Mat_MatTransMultDense *abt = c->abtdense; 1896faa55883SToby Isaac PetscErrorCode ierr; 1897faa55883SToby Isaac 1898faa55883SToby Isaac PetscFunctionBegin; 1899faa55883SToby Isaac switch (abt->alg) { 1900faa55883SToby Isaac case 1: 1901faa55883SToby Isaac ierr = MatMatTransposeMultNumeric_MPIDense_MPIDense_Cyclic(A, B, C);CHKERRQ(ierr); 1902faa55883SToby Isaac break; 1903faa55883SToby Isaac default: 1904faa55883SToby Isaac ierr = MatMatTransposeMultNumeric_MPIDense_MPIDense_Allgatherv(A, B, C);CHKERRQ(ierr); 1905faa55883SToby Isaac break; 1906faa55883SToby Isaac } 1907faa55883SToby Isaac PetscFunctionReturn(0); 1908faa55883SToby Isaac } 1909faa55883SToby Isaac 1910cc48ffa7SToby Isaac static PetscErrorCode MatMatTransposeMult_MPIDense_MPIDense(Mat A,Mat B, MatReuse scall, PetscReal fill, Mat *C) 1911cc48ffa7SToby Isaac { 1912cc48ffa7SToby Isaac PetscErrorCode ierr; 1913cc48ffa7SToby Isaac 1914cc48ffa7SToby Isaac PetscFunctionBegin; 1915cc48ffa7SToby Isaac if (scall == MAT_INITIAL_MATRIX) { 1916cc48ffa7SToby Isaac ierr = MatMatTransposeMultSymbolic_MPIDense_MPIDense(A,B,fill,C);CHKERRQ(ierr); 1917cc48ffa7SToby Isaac } 1918cc48ffa7SToby Isaac ierr = MatMatTransposeMultNumeric_MPIDense_MPIDense(A,B,*C);CHKERRQ(ierr); 1919cc48ffa7SToby Isaac PetscFunctionReturn(0); 1920cc48ffa7SToby Isaac } 1921cc48ffa7SToby Isaac 1922320f2790SHong Zhang PetscErrorCode MatDestroy_MatMatMult_MPIDense_MPIDense(Mat A) 1923320f2790SHong Zhang { 1924320f2790SHong Zhang PetscErrorCode ierr; 1925320f2790SHong Zhang Mat_MPIDense *a = (Mat_MPIDense*)A->data; 1926320f2790SHong Zhang Mat_MatMultDense *ab = a->abdense; 1927320f2790SHong Zhang 1928320f2790SHong Zhang PetscFunctionBegin; 1929320f2790SHong Zhang ierr = MatDestroy(&ab->Ce);CHKERRQ(ierr); 1930320f2790SHong Zhang ierr = MatDestroy(&ab->Ae);CHKERRQ(ierr); 1931320f2790SHong Zhang ierr = MatDestroy(&ab->Be);CHKERRQ(ierr); 1932320f2790SHong Zhang 1933320f2790SHong Zhang ierr = (ab->destroy)(A);CHKERRQ(ierr); 1934320f2790SHong Zhang ierr = PetscFree(ab);CHKERRQ(ierr); 1935320f2790SHong Zhang PetscFunctionReturn(0); 1936320f2790SHong Zhang } 1937320f2790SHong Zhang 19385242a7b1SHong Zhang #if defined(PETSC_HAVE_ELEMENTAL) 1939320f2790SHong Zhang PetscErrorCode MatMatMultNumeric_MPIDense_MPIDense(Mat A,Mat B,Mat C) 1940320f2790SHong Zhang { 1941320f2790SHong Zhang PetscErrorCode ierr; 1942320f2790SHong Zhang Mat_MPIDense *c=(Mat_MPIDense*)C->data; 1943320f2790SHong Zhang Mat_MatMultDense *ab=c->abdense; 1944320f2790SHong Zhang 1945320f2790SHong Zhang PetscFunctionBegin; 1946de0a22f0SHong Zhang ierr = MatConvert_MPIDense_Elemental(A,MATELEMENTAL,MAT_REUSE_MATRIX, &ab->Ae);CHKERRQ(ierr); 1947de0a22f0SHong Zhang ierr = MatConvert_MPIDense_Elemental(B,MATELEMENTAL,MAT_REUSE_MATRIX, &ab->Be);CHKERRQ(ierr); 1948320f2790SHong Zhang ierr = MatMatMultNumeric(ab->Ae,ab->Be,ab->Ce);CHKERRQ(ierr); 1949de0a22f0SHong Zhang ierr = MatConvert(ab->Ce,MATMPIDENSE,MAT_REUSE_MATRIX,&C);CHKERRQ(ierr); 1950320f2790SHong Zhang PetscFunctionReturn(0); 1951320f2790SHong Zhang } 1952320f2790SHong Zhang 1953320f2790SHong Zhang PetscErrorCode MatMatMultSymbolic_MPIDense_MPIDense(Mat A,Mat B,PetscReal fill,Mat *C) 1954320f2790SHong Zhang { 1955320f2790SHong Zhang PetscErrorCode ierr; 1956320f2790SHong Zhang Mat Ae,Be,Ce; 1957320f2790SHong Zhang Mat_MPIDense *c; 1958320f2790SHong Zhang Mat_MatMultDense *ab; 1959320f2790SHong Zhang 1960320f2790SHong Zhang PetscFunctionBegin; 1961320f2790SHong Zhang if (A->cmap->rstart != B->rmap->rstart || A->cmap->rend != B->rmap->rend) { 1962378336b6SHong 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); 1963320f2790SHong Zhang } 1964320f2790SHong Zhang 1965320f2790SHong Zhang /* convert A and B to Elemental matrices Ae and Be */ 1966320f2790SHong Zhang ierr = MatConvert(A,MATELEMENTAL,MAT_INITIAL_MATRIX, &Ae);CHKERRQ(ierr); 1967320f2790SHong Zhang ierr = MatConvert(B,MATELEMENTAL,MAT_INITIAL_MATRIX, &Be);CHKERRQ(ierr); 1968320f2790SHong Zhang 1969320f2790SHong Zhang /* Ce = Ae*Be */ 1970320f2790SHong Zhang ierr = MatMatMultSymbolic(Ae,Be,fill,&Ce);CHKERRQ(ierr); 1971320f2790SHong Zhang ierr = MatMatMultNumeric(Ae,Be,Ce);CHKERRQ(ierr); 1972320f2790SHong Zhang 1973320f2790SHong Zhang /* convert Ce to C */ 1974320f2790SHong Zhang ierr = MatConvert(Ce,MATMPIDENSE,MAT_INITIAL_MATRIX,C);CHKERRQ(ierr); 1975320f2790SHong Zhang 1976320f2790SHong Zhang /* create data structure for reuse Cdense */ 1977320f2790SHong Zhang ierr = PetscNew(&ab);CHKERRQ(ierr); 1978320f2790SHong Zhang c = (Mat_MPIDense*)(*C)->data; 1979320f2790SHong Zhang c->abdense = ab; 1980320f2790SHong Zhang 1981320f2790SHong Zhang ab->Ae = Ae; 1982320f2790SHong Zhang ab->Be = Be; 1983320f2790SHong Zhang ab->Ce = Ce; 1984320f2790SHong Zhang ab->destroy = (*C)->ops->destroy; 1985320f2790SHong Zhang (*C)->ops->destroy = MatDestroy_MatMatMult_MPIDense_MPIDense; 19869775878aSHong Zhang (*C)->ops->matmultnumeric = MatMatMultNumeric_MPIDense_MPIDense; 1987320f2790SHong Zhang PetscFunctionReturn(0); 1988320f2790SHong Zhang } 1989320f2790SHong Zhang 1990150d2497SBarry Smith PETSC_INTERN PetscErrorCode MatMatMult_MPIDense_MPIDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 1991320f2790SHong Zhang { 1992320f2790SHong Zhang PetscErrorCode ierr; 1993320f2790SHong Zhang 1994320f2790SHong Zhang PetscFunctionBegin; 199552c5f739Sprj- if (scall == MAT_INITIAL_MATRIX) { /* symbolic product includes numeric product */ 199657299f2fSHong Zhang ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr); 1997320f2790SHong Zhang ierr = MatMatMultSymbolic_MPIDense_MPIDense(A,B,fill,C);CHKERRQ(ierr); 199857299f2fSHong Zhang ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr); 1999320f2790SHong Zhang } else { 200057299f2fSHong Zhang ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr); 2001320f2790SHong Zhang ierr = MatMatMultNumeric_MPIDense_MPIDense(A,B,*C);CHKERRQ(ierr); 200257299f2fSHong Zhang ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr); 2003320f2790SHong Zhang } 2004320f2790SHong Zhang PetscFunctionReturn(0); 2005320f2790SHong Zhang } 20065242a7b1SHong Zhang #endif 200786aefd0dSHong Zhang 2008