1f4bdf6c4SBarry Smith 2f4bdf6c4SBarry Smith /* 3f4bdf6c4SBarry Smith Creates hypre ijmatrix from PETSc matrix 4f4bdf6c4SBarry Smith */ 5c6db04a5SJed Brown #include <petscsys.h> 607475bc1SBarry Smith #include <petsc-private/matimpl.h> 7c6db04a5SJed Brown #include <petscdmda.h> /*I "petscdmda.h" I*/ 8c6db04a5SJed Brown #include <../src/dm/impls/da/hypre/mhyp.h> 9f4bdf6c4SBarry Smith 10f4bdf6c4SBarry Smith #undef __FUNCT__ 11f4bdf6c4SBarry Smith #define __FUNCT__ "MatHYPRE_IJMatrixPreallocate" 12f4bdf6c4SBarry Smith PetscErrorCode MatHYPRE_IJMatrixPreallocate(Mat A_d, Mat A_o,HYPRE_IJMatrix ij) 13f4bdf6c4SBarry Smith { 14f4bdf6c4SBarry Smith PetscErrorCode ierr; 151a83f524SJed Brown PetscInt i,n_d,n_o; 161a83f524SJed Brown const PetscInt *ia_d,*ia_o; 17f4bdf6c4SBarry Smith PetscBool done_d=PETSC_FALSE,done_o=PETSC_FALSE; 18f4bdf6c4SBarry Smith PetscInt *nnz_d=PETSC_NULL,*nnz_o=PETSC_NULL; 19f4bdf6c4SBarry Smith 20f4bdf6c4SBarry Smith PetscFunctionBegin; 21f4bdf6c4SBarry Smith if (A_d) { /* determine number of nonzero entries in local diagonal part */ 22f4bdf6c4SBarry Smith ierr = MatGetRowIJ(A_d,0,PETSC_FALSE,PETSC_FALSE,&n_d,&ia_d,PETSC_NULL,&done_d);CHKERRQ(ierr); 23f4bdf6c4SBarry Smith if (done_d) { 24f4bdf6c4SBarry Smith ierr = PetscMalloc(n_d*sizeof(PetscInt),&nnz_d);CHKERRQ(ierr); 25f4bdf6c4SBarry Smith for (i=0; i<n_d; i++) { 26f4bdf6c4SBarry Smith nnz_d[i] = ia_d[i+1] - ia_d[i]; 27f4bdf6c4SBarry Smith } 28f4bdf6c4SBarry Smith } 29f4bdf6c4SBarry Smith ierr = MatRestoreRowIJ(A_d,0,PETSC_FALSE,PETSC_FALSE,&n_d,&ia_d,PETSC_NULL,&done_d);CHKERRQ(ierr); 30f4bdf6c4SBarry Smith } 31f4bdf6c4SBarry Smith if (A_o) { /* determine number of nonzero entries in local off-diagonal part */ 32f4bdf6c4SBarry Smith ierr = MatGetRowIJ(A_o,0,PETSC_FALSE,PETSC_FALSE,&n_o,&ia_o,PETSC_NULL,&done_o);CHKERRQ(ierr); 33f4bdf6c4SBarry Smith if (done_o) { 34f4bdf6c4SBarry Smith ierr = PetscMalloc(n_o*sizeof(PetscInt),&nnz_o);CHKERRQ(ierr); 35f4bdf6c4SBarry Smith for (i=0; i<n_o; i++) { 36f4bdf6c4SBarry Smith nnz_o[i] = ia_o[i+1] - ia_o[i]; 37f4bdf6c4SBarry Smith } 38f4bdf6c4SBarry Smith } 39f4bdf6c4SBarry Smith ierr = MatRestoreRowIJ(A_o,0,PETSC_FALSE,PETSC_FALSE,&n_o,&ia_o,PETSC_NULL,&done_o);CHKERRQ(ierr); 40f4bdf6c4SBarry Smith } 41f4bdf6c4SBarry Smith if (done_d) { /* set number of nonzeros in HYPRE IJ matrix */ 42f4bdf6c4SBarry Smith if (!done_o) { /* only diagonal part */ 43f4bdf6c4SBarry Smith ierr = PetscMalloc(n_d*sizeof(PetscInt),&nnz_o);CHKERRQ(ierr); 44f4bdf6c4SBarry Smith for (i=0; i<n_d; i++) { 45f4bdf6c4SBarry Smith nnz_o[i] = 0; 46f4bdf6c4SBarry Smith } 47f4bdf6c4SBarry Smith } 48fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_IJMatrixSetDiagOffdSizes,(ij,nnz_d,nnz_o)); 49f4bdf6c4SBarry Smith ierr = PetscFree(nnz_d);CHKERRQ(ierr); 50f4bdf6c4SBarry Smith ierr = PetscFree(nnz_o);CHKERRQ(ierr); 51f4bdf6c4SBarry Smith } 52f4bdf6c4SBarry Smith PetscFunctionReturn(0); 53f4bdf6c4SBarry Smith } 54f4bdf6c4SBarry Smith 55f4bdf6c4SBarry Smith #undef __FUNCT__ 56f4bdf6c4SBarry Smith #define __FUNCT__ "MatHYPRE_IJMatrixCreate" 57f4bdf6c4SBarry Smith PetscErrorCode MatHYPRE_IJMatrixCreate(Mat A,HYPRE_IJMatrix *ij) 58f4bdf6c4SBarry Smith { 59f4bdf6c4SBarry Smith PetscErrorCode ierr; 604ddd07fcSJed Brown PetscInt rstart,rend,cstart,cend; 61f4bdf6c4SBarry Smith 62f4bdf6c4SBarry Smith PetscFunctionBegin; 63f4bdf6c4SBarry Smith PetscValidHeaderSpecific(A,MAT_CLASSID,1); 64f4bdf6c4SBarry Smith PetscValidType(A,1); 65f4bdf6c4SBarry Smith PetscValidPointer(ij,2); 664994cf47SJed Brown ierr = MatSetUp(A);CHKERRQ(ierr); 67f4bdf6c4SBarry Smith rstart = A->rmap->rstart; 68f4bdf6c4SBarry Smith rend = A->rmap->rend; 69f4bdf6c4SBarry Smith cstart = A->cmap->rstart; 70f4bdf6c4SBarry Smith cend = A->cmap->rend; 71fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_IJMatrixCreate,(((PetscObject)A)->comm,rstart,rend-1,cstart,cend-1,ij)); 72fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_IJMatrixSetObjectType,(*ij,HYPRE_PARCSR)); 73f4bdf6c4SBarry Smith { 74f4bdf6c4SBarry Smith PetscBool same; 75f4bdf6c4SBarry Smith Mat A_d,A_o; 761a83f524SJed Brown const PetscInt *colmap; 77251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&same);CHKERRQ(ierr); 78f4bdf6c4SBarry Smith if (same) { 79f4bdf6c4SBarry Smith ierr = MatMPIAIJGetSeqAIJ(A,&A_d,&A_o,&colmap);CHKERRQ(ierr); 80f4bdf6c4SBarry Smith ierr = MatHYPRE_IJMatrixPreallocate(A_d,A_o,*ij);CHKERRQ(ierr); 81f4bdf6c4SBarry Smith PetscFunctionReturn(0); 82f4bdf6c4SBarry Smith } 83251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIBAIJ,&same);CHKERRQ(ierr); 84f4bdf6c4SBarry Smith if (same) { 85f4bdf6c4SBarry Smith ierr = MatMPIBAIJGetSeqBAIJ(A,&A_d,&A_o,&colmap);CHKERRQ(ierr); 86f4bdf6c4SBarry Smith ierr = MatHYPRE_IJMatrixPreallocate(A_d,A_o,*ij);CHKERRQ(ierr); 87f4bdf6c4SBarry Smith PetscFunctionReturn(0); 88f4bdf6c4SBarry Smith } 89251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&same);CHKERRQ(ierr); 90f4bdf6c4SBarry Smith if (same) { 91f4bdf6c4SBarry Smith ierr = MatHYPRE_IJMatrixPreallocate(A,PETSC_NULL,*ij);CHKERRQ(ierr); 92f4bdf6c4SBarry Smith PetscFunctionReturn(0); 93f4bdf6c4SBarry Smith } 94251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQBAIJ,&same);CHKERRQ(ierr); 95f4bdf6c4SBarry Smith if (same) { 96f4bdf6c4SBarry Smith ierr = MatHYPRE_IJMatrixPreallocate(A,PETSC_NULL,*ij);CHKERRQ(ierr); 97f4bdf6c4SBarry Smith PetscFunctionReturn(0); 98f4bdf6c4SBarry Smith } 99f4bdf6c4SBarry Smith } 100f4bdf6c4SBarry Smith PetscFunctionReturn(0); 101f4bdf6c4SBarry Smith } 102f4bdf6c4SBarry Smith 103f4bdf6c4SBarry Smith extern PetscErrorCode MatHYPRE_IJMatrixFastCopy_MPIAIJ(Mat,HYPRE_IJMatrix); 104f4bdf6c4SBarry Smith extern PetscErrorCode MatHYPRE_IJMatrixFastCopy_SeqAIJ(Mat,HYPRE_IJMatrix); 105f4bdf6c4SBarry Smith /* 106f4bdf6c4SBarry Smith Copies the data over (column indices, numerical values) to hypre matrix 107f4bdf6c4SBarry Smith */ 108f4bdf6c4SBarry Smith 109f4bdf6c4SBarry Smith #undef __FUNCT__ 110f4bdf6c4SBarry Smith #define __FUNCT__ "MatHYPRE_IJMatrixCopy" 111f4bdf6c4SBarry Smith PetscErrorCode MatHYPRE_IJMatrixCopy(Mat A,HYPRE_IJMatrix ij) 112f4bdf6c4SBarry Smith { 113f4bdf6c4SBarry Smith PetscErrorCode ierr; 114f4bdf6c4SBarry Smith PetscInt i,rstart,rend,ncols; 115f4bdf6c4SBarry Smith const PetscScalar *values; 116f4bdf6c4SBarry Smith const PetscInt *cols; 117f4bdf6c4SBarry Smith PetscBool flg; 118f4bdf6c4SBarry Smith 119f4bdf6c4SBarry Smith PetscFunctionBegin; 120251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&flg);CHKERRQ(ierr); 121f4bdf6c4SBarry Smith if (flg) { 122f4bdf6c4SBarry Smith ierr = MatHYPRE_IJMatrixFastCopy_MPIAIJ(A,ij);CHKERRQ(ierr); 123f4bdf6c4SBarry Smith PetscFunctionReturn(0); 124f4bdf6c4SBarry Smith } 125251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&flg);CHKERRQ(ierr); 126f4bdf6c4SBarry Smith if (flg) { 127f4bdf6c4SBarry Smith ierr = MatHYPRE_IJMatrixFastCopy_SeqAIJ(A,ij);CHKERRQ(ierr); 128f4bdf6c4SBarry Smith PetscFunctionReturn(0); 129f4bdf6c4SBarry Smith } 130f4bdf6c4SBarry Smith 131f4bdf6c4SBarry Smith ierr = PetscLogEventBegin(MAT_Convert,A,0,0,0);CHKERRQ(ierr); 132fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_IJMatrixInitialize,(ij)); 133f4bdf6c4SBarry Smith ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr); 134f4bdf6c4SBarry Smith for (i=rstart; i<rend; i++) { 135f4bdf6c4SBarry Smith ierr = MatGetRow(A,i,&ncols,&cols,&values);CHKERRQ(ierr); 136fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_IJMatrixSetValues,(ij,1,&ncols,&i,cols,values)); 137f4bdf6c4SBarry Smith ierr = MatRestoreRow(A,i,&ncols,&cols,&values);CHKERRQ(ierr); 138f4bdf6c4SBarry Smith } 139fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_IJMatrixAssemble,(ij)); 140f4bdf6c4SBarry Smith ierr = PetscLogEventEnd(MAT_Convert,A,0,0,0);CHKERRQ(ierr); 141f4bdf6c4SBarry Smith PetscFunctionReturn(0); 142f4bdf6c4SBarry Smith } 143f4bdf6c4SBarry Smith 144f4bdf6c4SBarry Smith /* 145f4bdf6c4SBarry Smith This copies the CSR format directly from the PETSc data structure to the hypre 146f4bdf6c4SBarry Smith data structure without calls to MatGetRow() or hypre's set values. 147f4bdf6c4SBarry Smith 148f4bdf6c4SBarry Smith */ 149c6db04a5SJed Brown #include <_hypre_IJ_mv.h> 150c6db04a5SJed Brown #include <HYPRE_IJ_mv.h> 151c6db04a5SJed Brown #include <../src/mat/impls/aij/mpi/mpiaij.h> 152f4bdf6c4SBarry Smith 153f4bdf6c4SBarry Smith #undef __FUNCT__ 154f4bdf6c4SBarry Smith #define __FUNCT__ "MatHYPRE_IJMatrixFastCopy_SeqAIJ" 155f4bdf6c4SBarry Smith PetscErrorCode MatHYPRE_IJMatrixFastCopy_SeqAIJ(Mat A,HYPRE_IJMatrix ij) 156f4bdf6c4SBarry Smith { 157f4bdf6c4SBarry Smith PetscErrorCode ierr; 158e497e08cSBarry Smith Mat_SeqAIJ *pdiag = (Mat_SeqAIJ*)A->data; 159f4bdf6c4SBarry Smith 160f4bdf6c4SBarry Smith hypre_ParCSRMatrix *par_matrix; 161f4bdf6c4SBarry Smith hypre_AuxParCSRMatrix *aux_matrix; 16225578ef6SJed Brown hypre_CSRMatrix *hdiag; 163f4bdf6c4SBarry Smith 164f4bdf6c4SBarry Smith PetscFunctionBegin; 165f4bdf6c4SBarry Smith PetscValidHeaderSpecific(A,MAT_CLASSID,1); 166f4bdf6c4SBarry Smith PetscValidType(A,1); 167f4bdf6c4SBarry Smith PetscValidPointer(ij,2); 168f4bdf6c4SBarry Smith 169f4bdf6c4SBarry Smith ierr = PetscLogEventBegin(MAT_Convert,A,0,0,0);CHKERRQ(ierr); 170fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_IJMatrixInitialize,(ij)); 171f4bdf6c4SBarry Smith par_matrix = (hypre_ParCSRMatrix*)hypre_IJMatrixObject(ij); 172f4bdf6c4SBarry Smith aux_matrix = (hypre_AuxParCSRMatrix*)hypre_IJMatrixTranslator(ij); 173f4bdf6c4SBarry Smith hdiag = hypre_ParCSRMatrixDiag(par_matrix); 174f4bdf6c4SBarry Smith 175f4bdf6c4SBarry Smith /* 176f4bdf6c4SBarry Smith this is the Hack part where we monkey directly with the hypre datastructures 177f4bdf6c4SBarry Smith */ 178f4bdf6c4SBarry Smith ierr = PetscMemcpy(hdiag->i,pdiag->i,(A->rmap->n + 1)*sizeof(PetscInt)); 179f4bdf6c4SBarry Smith ierr = PetscMemcpy(hdiag->j,pdiag->j,pdiag->nz*sizeof(PetscInt)); 180f4bdf6c4SBarry Smith ierr = PetscMemcpy(hdiag->data,pdiag->a,pdiag->nz*sizeof(PetscScalar)); 181f4bdf6c4SBarry Smith 182f4bdf6c4SBarry Smith hypre_AuxParCSRMatrixNeedAux(aux_matrix) = 0; 183fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_IJMatrixAssemble,(ij)); 184f4bdf6c4SBarry Smith ierr = PetscLogEventEnd(MAT_Convert,A,0,0,0);CHKERRQ(ierr); 185f4bdf6c4SBarry Smith PetscFunctionReturn(0); 186f4bdf6c4SBarry Smith } 187f4bdf6c4SBarry Smith 188f4bdf6c4SBarry Smith #undef __FUNCT__ 189f4bdf6c4SBarry Smith #define __FUNCT__ "MatHYPRE_IJMatrixFastCopy_MPIAIJ" 190f4bdf6c4SBarry Smith PetscErrorCode MatHYPRE_IJMatrixFastCopy_MPIAIJ(Mat A,HYPRE_IJMatrix ij) 191f4bdf6c4SBarry Smith { 192f4bdf6c4SBarry Smith PetscErrorCode ierr; 193f4bdf6c4SBarry Smith Mat_MPIAIJ *pA = (Mat_MPIAIJ*)A->data; 194f4bdf6c4SBarry Smith Mat_SeqAIJ *pdiag,*poffd; 195f4bdf6c4SBarry Smith PetscInt i,*garray = pA->garray,*jj,cstart,*pjj; 196f4bdf6c4SBarry Smith 197f4bdf6c4SBarry Smith hypre_ParCSRMatrix *par_matrix; 198f4bdf6c4SBarry Smith hypre_AuxParCSRMatrix *aux_matrix; 199f4bdf6c4SBarry Smith hypre_CSRMatrix *hdiag,*hoffd; 200f4bdf6c4SBarry Smith 201f4bdf6c4SBarry Smith PetscFunctionBegin; 202f4bdf6c4SBarry Smith PetscValidHeaderSpecific(A,MAT_CLASSID,1); 203f4bdf6c4SBarry Smith PetscValidType(A,1); 204f4bdf6c4SBarry Smith PetscValidPointer(ij,2); 205f4bdf6c4SBarry Smith pdiag = (Mat_SeqAIJ*) pA->A->data; 206f4bdf6c4SBarry Smith poffd = (Mat_SeqAIJ*) pA->B->data; 207f4bdf6c4SBarry Smith /* cstart is only valid for square MPIAIJ layed out in the usual way */ 208f4bdf6c4SBarry Smith ierr = MatGetOwnershipRange(A,&cstart,PETSC_NULL);CHKERRQ(ierr); 209f4bdf6c4SBarry Smith 210f4bdf6c4SBarry Smith ierr = PetscLogEventBegin(MAT_Convert,A,0,0,0);CHKERRQ(ierr); 211fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_IJMatrixInitialize,(ij)); 212f4bdf6c4SBarry Smith par_matrix = (hypre_ParCSRMatrix*)hypre_IJMatrixObject(ij); 213f4bdf6c4SBarry Smith aux_matrix = (hypre_AuxParCSRMatrix*)hypre_IJMatrixTranslator(ij); 214f4bdf6c4SBarry Smith hdiag = hypre_ParCSRMatrixDiag(par_matrix); 215f4bdf6c4SBarry Smith hoffd = hypre_ParCSRMatrixOffd(par_matrix); 216f4bdf6c4SBarry Smith 217f4bdf6c4SBarry Smith /* 218f4bdf6c4SBarry Smith this is the Hack part where we monkey directly with the hypre datastructures 219f4bdf6c4SBarry Smith */ 220f4bdf6c4SBarry Smith ierr = PetscMemcpy(hdiag->i,pdiag->i,(pA->A->rmap->n + 1)*sizeof(PetscInt)); 221f4bdf6c4SBarry Smith /* need to shift the diag column indices (hdiag->j) back to global numbering since hypre is expecting this */ 222f4bdf6c4SBarry Smith jj = hdiag->j; 223f4bdf6c4SBarry Smith pjj = pdiag->j; 224*8865f1eaSKarl Rupp for (i=0; i<pdiag->nz; i++) jj[i] = cstart + pjj[i]; 225f4bdf6c4SBarry Smith ierr = PetscMemcpy(hdiag->data,pdiag->a,pdiag->nz*sizeof(PetscScalar)); 226f4bdf6c4SBarry Smith ierr = PetscMemcpy(hoffd->i,poffd->i,(pA->A->rmap->n + 1)*sizeof(PetscInt)); 227f4bdf6c4SBarry Smith /* need to move the offd column indices (hoffd->j) back to global numbering since hypre is expecting this 228f4bdf6c4SBarry Smith If we hacked a hypre a bit more we might be able to avoid this step */ 229f4bdf6c4SBarry Smith jj = hoffd->j; 230f4bdf6c4SBarry Smith pjj = poffd->j; 231*8865f1eaSKarl Rupp for (i=0; i<poffd->nz; i++) jj[i] = garray[pjj[i]]; 232f4bdf6c4SBarry Smith ierr = PetscMemcpy(hoffd->data,poffd->a,poffd->nz*sizeof(PetscScalar)); 233f4bdf6c4SBarry Smith 234f4bdf6c4SBarry Smith hypre_AuxParCSRMatrixNeedAux(aux_matrix) = 0; 235fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_IJMatrixAssemble,(ij)); 236f4bdf6c4SBarry Smith ierr = PetscLogEventEnd(MAT_Convert,A,0,0,0);CHKERRQ(ierr); 237f4bdf6c4SBarry Smith PetscFunctionReturn(0); 238f4bdf6c4SBarry Smith } 239f4bdf6c4SBarry Smith 240f4bdf6c4SBarry Smith /* 241f4bdf6c4SBarry Smith Does NOT copy the data over, instead uses DIRECTLY the pointers from the PETSc MPIAIJ format 242f4bdf6c4SBarry Smith 243f4bdf6c4SBarry Smith This is UNFINISHED and does NOT work! The problem is that hypre puts the diagonal entry first 244f4bdf6c4SBarry Smith which will corrupt the PETSc data structure if we did this. Need a work around to this problem. 245f4bdf6c4SBarry Smith */ 246c6db04a5SJed Brown #include <_hypre_IJ_mv.h> 247c6db04a5SJed Brown #include <HYPRE_IJ_mv.h> 248f4bdf6c4SBarry Smith 249f4bdf6c4SBarry Smith #undef __FUNCT__ 250f4bdf6c4SBarry Smith #define __FUNCT__ "MatHYPRE_IJMatrixLink" 251f4bdf6c4SBarry Smith PetscErrorCode MatHYPRE_IJMatrixLink(Mat A,HYPRE_IJMatrix *ij) 252f4bdf6c4SBarry Smith { 253f4bdf6c4SBarry Smith PetscErrorCode ierr; 2544ddd07fcSJed Brown PetscInt rstart,rend,cstart,cend; 255f4bdf6c4SBarry Smith PetscBool flg; 256f4bdf6c4SBarry Smith hypre_AuxParCSRMatrix *aux_matrix; 257f4bdf6c4SBarry Smith 258f4bdf6c4SBarry Smith PetscFunctionBegin; 259f4bdf6c4SBarry Smith PetscValidHeaderSpecific(A,MAT_CLASSID,1); 260f4bdf6c4SBarry Smith PetscValidType(A,1); 261f4bdf6c4SBarry Smith PetscValidPointer(ij,2); 262251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&flg);CHKERRQ(ierr); 263f4bdf6c4SBarry Smith if (!flg) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"Can only use with PETSc MPIAIJ matrices"); 2644994cf47SJed Brown ierr = MatSetUp(A);CHKERRQ(ierr); 265f4bdf6c4SBarry Smith 266f4bdf6c4SBarry Smith ierr = PetscLogEventBegin(MAT_Convert,A,0,0,0);CHKERRQ(ierr); 267f4bdf6c4SBarry Smith rstart = A->rmap->rstart; 268f4bdf6c4SBarry Smith rend = A->rmap->rend; 269f4bdf6c4SBarry Smith cstart = A->cmap->rstart; 270f4bdf6c4SBarry Smith cend = A->cmap->rend; 271fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_IJMatrixCreate,(((PetscObject)A)->comm,rstart,rend-1,cstart,cend-1,ij)); 272fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_IJMatrixSetObjectType,(*ij,HYPRE_PARCSR)); 273f4bdf6c4SBarry Smith 274fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_IJMatrixInitialize,(*ij)); 275fd3f9acdSBarry Smith PetscStackCall("hypre_IJMatrixTranslator",aux_matrix = (hypre_AuxParCSRMatrix*)hypre_IJMatrixTranslator(*ij)); 276f4bdf6c4SBarry Smith 277f4bdf6c4SBarry Smith hypre_AuxParCSRMatrixNeedAux(aux_matrix) = 0; 278f4bdf6c4SBarry Smith 279f4bdf6c4SBarry Smith /* this is the Hack part where we monkey directly with the hypre datastructures */ 280f4bdf6c4SBarry Smith 281fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_IJMatrixAssemble,(*ij)); 282f4bdf6c4SBarry Smith ierr = PetscLogEventEnd(MAT_Convert,A,0,0,0);CHKERRQ(ierr); 283f4bdf6c4SBarry Smith PetscFunctionReturn(0); 284f4bdf6c4SBarry Smith } 285f4bdf6c4SBarry Smith 286f4bdf6c4SBarry Smith /* -----------------------------------------------------------------------------------------------------------------*/ 287f4bdf6c4SBarry Smith 288f4bdf6c4SBarry Smith /*MC 289f4bdf6c4SBarry Smith MATHYPRESTRUCT - MATHYPRESTRUCT = "hyprestruct" - A matrix type to be used for parallel sparse matrices 290f4bdf6c4SBarry Smith based on the hypre HYPRE_StructMatrix. 291f4bdf6c4SBarry Smith 292f4bdf6c4SBarry Smith Level: intermediate 293f4bdf6c4SBarry Smith 294f4bdf6c4SBarry Smith Notes: Unlike the more general support for blocks in hypre this allows only one block per process and requires the block 295f4bdf6c4SBarry Smith be defined by a DMDA. 296f4bdf6c4SBarry Smith 297c688c046SMatthew G Knepley The matrix needs a DMDA associated with it by either a call to MatSetupDM() or if the matrix is obtained from DMCreateMatrix() 298f4bdf6c4SBarry Smith 299c688c046SMatthew G Knepley .seealso: MatCreate(), PCPFMG, MatSetupDM(), DMCreateMatrix() 300f4bdf6c4SBarry Smith M*/ 301f4bdf6c4SBarry Smith 302f4bdf6c4SBarry Smith 303f4bdf6c4SBarry Smith #undef __FUNCT__ 304f4bdf6c4SBarry Smith #define __FUNCT__ "MatSetValuesLocal_HYPREStruct_3d" 305f4bdf6c4SBarry Smith PetscErrorCode MatSetValuesLocal_HYPREStruct_3d(Mat mat,PetscInt nrow,const PetscInt irow[],PetscInt ncol,const PetscInt icol[],const PetscScalar y[],InsertMode addv) 306f4bdf6c4SBarry Smith { 307f4bdf6c4SBarry Smith PetscErrorCode ierr; 308f4bdf6c4SBarry Smith PetscInt i,j,stencil,index[3],row,entries[7]; 309f4bdf6c4SBarry Smith const PetscScalar *values = y; 310f4bdf6c4SBarry Smith Mat_HYPREStruct *ex = (Mat_HYPREStruct*) mat->data; 311f4bdf6c4SBarry Smith 312f4bdf6c4SBarry Smith PetscFunctionBegin; 313f4bdf6c4SBarry Smith for (i=0; i<nrow; i++) { 314f4bdf6c4SBarry Smith for (j=0; j<ncol; j++) { 315f4bdf6c4SBarry Smith stencil = icol[j] - irow[i]; 316f4bdf6c4SBarry Smith if (!stencil) { 317f4bdf6c4SBarry Smith entries[j] = 3; 318f4bdf6c4SBarry Smith } else if (stencil == -1) { 319f4bdf6c4SBarry Smith entries[j] = 2; 320f4bdf6c4SBarry Smith } else if (stencil == 1) { 321f4bdf6c4SBarry Smith entries[j] = 4; 322f4bdf6c4SBarry Smith } else if (stencil == -ex->gnx) { 323f4bdf6c4SBarry Smith entries[j] = 1; 324f4bdf6c4SBarry Smith } else if (stencil == ex->gnx) { 325f4bdf6c4SBarry Smith entries[j] = 5; 326f4bdf6c4SBarry Smith } else if (stencil == -ex->gnxgny) { 327f4bdf6c4SBarry Smith entries[j] = 0; 328f4bdf6c4SBarry Smith } else if (stencil == ex->gnxgny) { 329f4bdf6c4SBarry Smith entries[j] = 6; 330f4bdf6c4SBarry Smith } else SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Local row %D local column %D have bad stencil %D",irow[i],icol[j],stencil); 331f4bdf6c4SBarry Smith } 332f4bdf6c4SBarry Smith row = ex->gindices[irow[i]] - ex->rstart; 333f4bdf6c4SBarry Smith index[0] = ex->xs + (row % ex->nx); 334f4bdf6c4SBarry Smith index[1] = ex->ys + ((row/ex->nx) % ex->ny); 335f4bdf6c4SBarry Smith index[2] = ex->zs + (row/(ex->nxny)); 336f4bdf6c4SBarry Smith if (addv == ADD_VALUES) { 337fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_StructMatrixAddToValues,(ex->hmat,index,ncol,entries,(PetscScalar*)values)); 338f4bdf6c4SBarry Smith } else { 339fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_StructMatrixSetValues,(ex->hmat,index,ncol,entries,(PetscScalar*)values)); 340f4bdf6c4SBarry Smith } 341f4bdf6c4SBarry Smith values += ncol; 342f4bdf6c4SBarry Smith } 343f4bdf6c4SBarry Smith PetscFunctionReturn(0); 344f4bdf6c4SBarry Smith } 345f4bdf6c4SBarry Smith 346f4bdf6c4SBarry Smith #undef __FUNCT__ 347f4bdf6c4SBarry Smith #define __FUNCT__ "MatZeroRowsLocal_HYPREStruct_3d" 348f4bdf6c4SBarry Smith PetscErrorCode MatZeroRowsLocal_HYPREStruct_3d(Mat mat,PetscInt nrow,const PetscInt irow[],PetscScalar d,Vec x,Vec b) 349f4bdf6c4SBarry Smith { 350f4bdf6c4SBarry Smith PetscErrorCode ierr; 351f4bdf6c4SBarry Smith PetscInt i,index[3],row,entries[7] = {0,1,2,3,4,5,6}; 352f4bdf6c4SBarry Smith PetscScalar values[7]; 353f4bdf6c4SBarry Smith Mat_HYPREStruct *ex = (Mat_HYPREStruct*) mat->data; 354f4bdf6c4SBarry Smith 355f4bdf6c4SBarry Smith PetscFunctionBegin; 356f4bdf6c4SBarry Smith if (x && b) SETERRQ(((PetscObject)mat)->comm,PETSC_ERR_SUP,"No support"); 357f4bdf6c4SBarry Smith ierr = PetscMemzero(values,7*sizeof(PetscScalar));CHKERRQ(ierr); 358f4bdf6c4SBarry Smith values[3] = d; 359f4bdf6c4SBarry Smith for (i=0; i<nrow; i++) { 360f4bdf6c4SBarry Smith row = ex->gindices[irow[i]] - ex->rstart; 361f4bdf6c4SBarry Smith index[0] = ex->xs + (row % ex->nx); 362f4bdf6c4SBarry Smith index[1] = ex->ys + ((row/ex->nx) % ex->ny); 363f4bdf6c4SBarry Smith index[2] = ex->zs + (row/(ex->nxny)); 364fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_StructMatrixSetValues,(ex->hmat,index,7,entries,values)); 365f4bdf6c4SBarry Smith } 366fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_StructMatrixAssemble,(ex->hmat)); 367f4bdf6c4SBarry Smith PetscFunctionReturn(0); 368f4bdf6c4SBarry Smith } 369f4bdf6c4SBarry Smith 370f4bdf6c4SBarry Smith #undef __FUNCT__ 371f4bdf6c4SBarry Smith #define __FUNCT__ "MatZeroEntries_HYPREStruct_3d" 372f4bdf6c4SBarry Smith PetscErrorCode MatZeroEntries_HYPREStruct_3d(Mat mat) 373f4bdf6c4SBarry Smith { 374f4bdf6c4SBarry Smith PetscErrorCode ierr; 375f4bdf6c4SBarry Smith PetscInt indices[7] = {0,1,2,3,4,5,6}; 376f4bdf6c4SBarry Smith Mat_HYPREStruct *ex = (Mat_HYPREStruct*) mat->data; 377f4bdf6c4SBarry Smith 378f4bdf6c4SBarry Smith PetscFunctionBegin; 379f4bdf6c4SBarry Smith /* hypre has no public interface to do this */ 380fd3f9acdSBarry Smith PetscStackCallStandard(hypre_StructMatrixClearBoxValues,(ex->hmat,&ex->hbox,7,indices,0,1)); 381fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_StructMatrixAssemble,(ex->hmat)); 382f4bdf6c4SBarry Smith PetscFunctionReturn(0); 383f4bdf6c4SBarry Smith } 384f4bdf6c4SBarry Smith 385f4bdf6c4SBarry Smith #undef __FUNCT__ 386c688c046SMatthew G Knepley #define __FUNCT__ "MatSetupDM_HYPREStruct" 387c688c046SMatthew G Knepley PetscErrorCode MatSetupDM_HYPREStruct(Mat mat,DM da) 388f4bdf6c4SBarry Smith { 389f4bdf6c4SBarry Smith PetscErrorCode ierr; 390f4bdf6c4SBarry Smith Mat_HYPREStruct *ex = (Mat_HYPREStruct*) mat->data; 3914ddd07fcSJed Brown PetscInt dim,dof,sw[3],nx,ny,nz,ilower[3],iupper[3],ssize,i; 3927ab44359SJed Brown DMDABoundaryType px,py,pz; 393f4bdf6c4SBarry Smith DMDAStencilType st; 394f4bdf6c4SBarry Smith 395f4bdf6c4SBarry Smith PetscFunctionBegin; 396f4bdf6c4SBarry Smith ex->da = da; 397f4bdf6c4SBarry Smith ierr = PetscObjectReference((PetscObject)da);CHKERRQ(ierr); 398f4bdf6c4SBarry Smith 3997ab44359SJed Brown ierr = DMDAGetInfo(ex->da,&dim,0,0,0,0,0,0,&dof,&sw[0],&px,&py,&pz,&st);CHKERRQ(ierr); 400f4bdf6c4SBarry Smith ierr = DMDAGetCorners(ex->da,&ilower[0],&ilower[1],&ilower[2],&iupper[0],&iupper[1],&iupper[2]);CHKERRQ(ierr); 401f4bdf6c4SBarry Smith iupper[0] += ilower[0] - 1; 402f4bdf6c4SBarry Smith iupper[1] += ilower[1] - 1; 403f4bdf6c4SBarry Smith iupper[2] += ilower[2] - 1; 404f4bdf6c4SBarry Smith 405f4bdf6c4SBarry Smith /* the hypre_Box is used to zero out the matrix entries in MatZeroValues() */ 406f4bdf6c4SBarry Smith ex->hbox.imin[0] = ilower[0]; 407f4bdf6c4SBarry Smith ex->hbox.imin[1] = ilower[1]; 408f4bdf6c4SBarry Smith ex->hbox.imin[2] = ilower[2]; 409f4bdf6c4SBarry Smith ex->hbox.imax[0] = iupper[0]; 410f4bdf6c4SBarry Smith ex->hbox.imax[1] = iupper[1]; 411f4bdf6c4SBarry Smith ex->hbox.imax[2] = iupper[2]; 412f4bdf6c4SBarry Smith 413f4bdf6c4SBarry Smith /* create the hypre grid object and set its information */ 414f4bdf6c4SBarry Smith if (dof > 1) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"Currently only support for scalar problems"); 4157ab44359SJed Brown if (px || py || pz) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"Ask us to add periodic support by calling HYPRE_StructGridSetPeriodic()"); 416fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_StructGridCreate,(ex->hcomm,dim,&ex->hgrid)); 417fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_StructGridSetExtents,(ex->hgrid,ilower,iupper)); 418fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_StructGridAssemble,(ex->hgrid)); 419f4bdf6c4SBarry Smith 420f4bdf6c4SBarry Smith sw[1] = sw[0]; 421f4bdf6c4SBarry Smith sw[2] = sw[1]; 422fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_StructGridSetNumGhost,(ex->hgrid,sw)); 423f4bdf6c4SBarry Smith 424f4bdf6c4SBarry Smith /* create the hypre stencil object and set its information */ 425f4bdf6c4SBarry Smith if (sw[0] > 1) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"Ask us to add support for wider stencils"); 426f4bdf6c4SBarry Smith if (st == DMDA_STENCIL_BOX) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"Ask us to add support for box stencils"); 427f4bdf6c4SBarry Smith if (dim == 1) { 4284ddd07fcSJed Brown PetscInt offsets[3][1] = {{-1},{0},{1}}; 429f4bdf6c4SBarry Smith ssize = 3; 430fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_StructStencilCreate,(dim,ssize,&ex->hstencil)); 431f4bdf6c4SBarry Smith for (i=0; i<ssize; i++) { 432fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_StructStencilSetElement,(ex->hstencil,i,offsets[i])); 433f4bdf6c4SBarry Smith } 434f4bdf6c4SBarry Smith } else if (dim == 2) { 4354ddd07fcSJed Brown PetscInt offsets[5][2] = {{0,-1},{-1,0},{0,0},{1,0},{0,1}}; 436f4bdf6c4SBarry Smith ssize = 5; 437fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_StructStencilCreate,(dim,ssize,&ex->hstencil)); 438f4bdf6c4SBarry Smith for (i=0; i<ssize; i++) { 439fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_StructStencilSetElement,(ex->hstencil,i,offsets[i])); 440f4bdf6c4SBarry Smith } 441f4bdf6c4SBarry Smith } else if (dim == 3) { 4424ddd07fcSJed Brown PetscInt offsets[7][3] = {{0,0,-1},{0,-1,0},{-1,0,0},{0,0,0},{1,0,0},{0,1,0},{0,0,1}}; 443f4bdf6c4SBarry Smith ssize = 7; 444fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_StructStencilCreate,(dim,ssize,&ex->hstencil)); 445f4bdf6c4SBarry Smith for (i=0; i<ssize; i++) { 446fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_StructStencilSetElement,(ex->hstencil,i,offsets[i])); 447f4bdf6c4SBarry Smith } 448f4bdf6c4SBarry Smith } 449f4bdf6c4SBarry Smith 450f4bdf6c4SBarry Smith /* create the HYPRE vector for rhs and solution */ 451fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_StructVectorCreate,(ex->hcomm,ex->hgrid,&ex->hb)); 452fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_StructVectorCreate,(ex->hcomm,ex->hgrid,&ex->hx)); 453fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_StructVectorInitialize,(ex->hb)); 454fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_StructVectorInitialize,(ex->hx)); 455fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_StructVectorAssemble,(ex->hb)); 456fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_StructVectorAssemble,(ex->hx)); 457f4bdf6c4SBarry Smith 458f4bdf6c4SBarry Smith /* create the hypre matrix object and set its information */ 459fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_StructMatrixCreate,(ex->hcomm,ex->hgrid,ex->hstencil,&ex->hmat)); 460fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_StructGridDestroy,(ex->hgrid)); 461fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_StructStencilDestroy,(ex->hstencil)); 462f4bdf6c4SBarry Smith if (ex->needsinitialization) { 463fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_StructMatrixInitialize,(ex->hmat)); 464f4bdf6c4SBarry Smith ex->needsinitialization = PETSC_FALSE; 465f4bdf6c4SBarry Smith } 466f4bdf6c4SBarry Smith 467f4bdf6c4SBarry Smith /* set the global and local sizes of the matrix */ 468f4bdf6c4SBarry Smith ierr = DMDAGetCorners(da,0,0,0,&nx,&ny,&nz);CHKERRQ(ierr); 469f4bdf6c4SBarry Smith ierr = MatSetSizes(mat,dof*nx*ny*nz,dof*nx*ny*nz,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr); 470f4bdf6c4SBarry Smith ierr = PetscLayoutSetBlockSize(mat->rmap,1);CHKERRQ(ierr); 471f4bdf6c4SBarry Smith ierr = PetscLayoutSetBlockSize(mat->cmap,1);CHKERRQ(ierr); 472f4bdf6c4SBarry Smith ierr = PetscLayoutSetUp(mat->rmap);CHKERRQ(ierr); 473f4bdf6c4SBarry Smith ierr = PetscLayoutSetUp(mat->cmap);CHKERRQ(ierr); 474f4bdf6c4SBarry Smith 475f4bdf6c4SBarry Smith if (dim == 3) { 476f4bdf6c4SBarry Smith mat->ops->setvalueslocal = MatSetValuesLocal_HYPREStruct_3d; 477f4bdf6c4SBarry Smith mat->ops->zerorowslocal = MatZeroRowsLocal_HYPREStruct_3d; 478f4bdf6c4SBarry Smith mat->ops->zeroentries = MatZeroEntries_HYPREStruct_3d; 479*8865f1eaSKarl Rupp 480f4bdf6c4SBarry Smith ierr = MatZeroEntries_HYPREStruct_3d(mat);CHKERRQ(ierr); 481f4bdf6c4SBarry Smith } else SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"Only support for 3d DMDA currently"); 482f4bdf6c4SBarry Smith 483f4bdf6c4SBarry Smith /* get values that will be used repeatedly in MatSetValuesLocal() and MatZeroRowsLocal() repeatedly */ 484f4bdf6c4SBarry Smith ierr = MatGetOwnershipRange(mat,&ex->rstart,PETSC_NULL);CHKERRQ(ierr); 485f4bdf6c4SBarry Smith ierr = DMDAGetGlobalIndices(ex->da,PETSC_NULL,&ex->gindices);CHKERRQ(ierr); 486f4bdf6c4SBarry Smith ierr = DMDAGetGhostCorners(ex->da,0,0,0,&ex->gnx,&ex->gnxgny,0);CHKERRQ(ierr); 487f4bdf6c4SBarry Smith ex->gnxgny *= ex->gnx; 488f4bdf6c4SBarry Smith ierr = DMDAGetCorners(ex->da,&ex->xs,&ex->ys,&ex->zs,&ex->nx,&ex->ny,0);CHKERRQ(ierr); 489f4bdf6c4SBarry Smith ex->nxny = ex->nx*ex->ny; 490f4bdf6c4SBarry Smith PetscFunctionReturn(0); 491f4bdf6c4SBarry Smith } 492f4bdf6c4SBarry Smith 493f4bdf6c4SBarry Smith #undef __FUNCT__ 494f4bdf6c4SBarry Smith #define __FUNCT__ "MatMult_HYPREStruct" 495f4bdf6c4SBarry Smith PetscErrorCode MatMult_HYPREStruct(Mat A,Vec x,Vec y) 496f4bdf6c4SBarry Smith { 497f4bdf6c4SBarry Smith PetscErrorCode ierr; 498f4bdf6c4SBarry Smith PetscScalar *xx,*yy; 4994ddd07fcSJed Brown PetscInt ilower[3],iupper[3]; 500f4bdf6c4SBarry Smith Mat_HYPREStruct *mx = (Mat_HYPREStruct*)(A->data); 501f4bdf6c4SBarry Smith 502f4bdf6c4SBarry Smith PetscFunctionBegin; 503f4bdf6c4SBarry Smith ierr = DMDAGetCorners(mx->da,&ilower[0],&ilower[1],&ilower[2],&iupper[0],&iupper[1],&iupper[2]);CHKERRQ(ierr); 504*8865f1eaSKarl Rupp 505f4bdf6c4SBarry Smith iupper[0] += ilower[0] - 1; 506f4bdf6c4SBarry Smith iupper[1] += ilower[1] - 1; 507f4bdf6c4SBarry Smith iupper[2] += ilower[2] - 1; 508f4bdf6c4SBarry Smith 509f4bdf6c4SBarry Smith /* copy x values over to hypre */ 510fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_StructVectorSetConstantValues,(mx->hb,0.0)); 511f4bdf6c4SBarry Smith ierr = VecGetArray(x,&xx);CHKERRQ(ierr); 512fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_StructVectorSetBoxValues,(mx->hb,ilower,iupper,xx)); 513f4bdf6c4SBarry Smith ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr); 514fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_StructVectorAssemble,(mx->hb)); 515fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_StructMatrixMatvec,(1.0,mx->hmat,mx->hb,0.0,mx->hx)); 516f4bdf6c4SBarry Smith 517f4bdf6c4SBarry Smith /* copy solution values back to PETSc */ 518f4bdf6c4SBarry Smith ierr = VecGetArray(y,&yy);CHKERRQ(ierr); 519fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_StructVectorGetBoxValues,(mx->hx,ilower,iupper,yy)); 520f4bdf6c4SBarry Smith ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr); 521f4bdf6c4SBarry Smith PetscFunctionReturn(0); 522f4bdf6c4SBarry Smith } 523f4bdf6c4SBarry Smith 524f4bdf6c4SBarry Smith #undef __FUNCT__ 525f4bdf6c4SBarry Smith #define __FUNCT__ "MatAssemblyEnd_HYPREStruct" 526f4bdf6c4SBarry Smith PetscErrorCode MatAssemblyEnd_HYPREStruct(Mat mat,MatAssemblyType mode) 527f4bdf6c4SBarry Smith { 528f4bdf6c4SBarry Smith Mat_HYPREStruct *ex = (Mat_HYPREStruct*) mat->data; 529f4bdf6c4SBarry Smith PetscErrorCode ierr; 530f4bdf6c4SBarry Smith 531f4bdf6c4SBarry Smith PetscFunctionBegin; 532fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_StructMatrixAssemble,(ex->hmat)); 533fd3f9acdSBarry Smith /* PetscStackCallStandard(HYPRE_StructMatrixPrint,("dummy",ex->hmat,0)); */ 534f4bdf6c4SBarry Smith PetscFunctionReturn(0); 535f4bdf6c4SBarry Smith } 536f4bdf6c4SBarry Smith 537f4bdf6c4SBarry Smith #undef __FUNCT__ 538f4bdf6c4SBarry Smith #define __FUNCT__ "MatZeroEntries_HYPREStruct" 539f4bdf6c4SBarry Smith PetscErrorCode MatZeroEntries_HYPREStruct(Mat mat) 540f4bdf6c4SBarry Smith { 541f4bdf6c4SBarry Smith PetscFunctionBegin; 542f4bdf6c4SBarry Smith /* before the DMDA is set to the matrix the zero doesn't need to do anything */ 543f4bdf6c4SBarry Smith PetscFunctionReturn(0); 544f4bdf6c4SBarry Smith } 545f4bdf6c4SBarry Smith 546f4bdf6c4SBarry Smith 547f4bdf6c4SBarry Smith #undef __FUNCT__ 548f4bdf6c4SBarry Smith #define __FUNCT__ "MatDestroy_HYPREStruct" 549f4bdf6c4SBarry Smith PetscErrorCode MatDestroy_HYPREStruct(Mat mat) 550f4bdf6c4SBarry Smith { 551f4bdf6c4SBarry Smith Mat_HYPREStruct *ex = (Mat_HYPREStruct*) mat->data; 552f4bdf6c4SBarry Smith PetscErrorCode ierr; 553f4bdf6c4SBarry Smith 554f4bdf6c4SBarry Smith PetscFunctionBegin; 555fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_StructMatrixDestroy,(ex->hmat)); 556fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_StructVectorDestroy,(ex->hx)); 557fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_StructVectorDestroy,(ex->hb)); 558f4bdf6c4SBarry Smith PetscFunctionReturn(0); 559f4bdf6c4SBarry Smith } 560f4bdf6c4SBarry Smith 561f4bdf6c4SBarry Smith 562f4bdf6c4SBarry Smith EXTERN_C_BEGIN 563f4bdf6c4SBarry Smith #undef __FUNCT__ 564f4bdf6c4SBarry Smith #define __FUNCT__ "MatCreate_HYPREStruct" 565f4bdf6c4SBarry Smith PetscErrorCode MatCreate_HYPREStruct(Mat B) 566f4bdf6c4SBarry Smith { 567f4bdf6c4SBarry Smith Mat_HYPREStruct *ex; 568f4bdf6c4SBarry Smith PetscErrorCode ierr; 569f4bdf6c4SBarry Smith 570f4bdf6c4SBarry Smith PetscFunctionBegin; 571f4bdf6c4SBarry Smith ierr = PetscNewLog(B,Mat_HYPREStruct,&ex);CHKERRQ(ierr); 572f4bdf6c4SBarry Smith B->data = (void*)ex; 573f4bdf6c4SBarry Smith B->rmap->bs = 1; 574f4bdf6c4SBarry Smith B->assembled = PETSC_FALSE; 575f4bdf6c4SBarry Smith 576f4bdf6c4SBarry Smith B->insertmode = NOT_SET_VALUES; 577f4bdf6c4SBarry Smith 578f4bdf6c4SBarry Smith B->ops->assemblyend = MatAssemblyEnd_HYPREStruct; 579f4bdf6c4SBarry Smith B->ops->mult = MatMult_HYPREStruct; 580f4bdf6c4SBarry Smith B->ops->zeroentries = MatZeroEntries_HYPREStruct; 581f4bdf6c4SBarry Smith B->ops->destroy = MatDestroy_HYPREStruct; 582f4bdf6c4SBarry Smith 583f4bdf6c4SBarry Smith ex->needsinitialization = PETSC_TRUE; 584f4bdf6c4SBarry Smith 585f4bdf6c4SBarry Smith ierr = MPI_Comm_dup(((PetscObject)B)->comm,&(ex->hcomm));CHKERRQ(ierr); 586c688c046SMatthew G Knepley ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSetupDM_C","MatSetupDM_HYPREStruct",MatSetupDM_HYPREStruct);CHKERRQ(ierr); 587f4bdf6c4SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)B,MATHYPRESTRUCT);CHKERRQ(ierr); 588f4bdf6c4SBarry Smith PetscFunctionReturn(0); 589f4bdf6c4SBarry Smith } 590f4bdf6c4SBarry Smith EXTERN_C_END 591f4bdf6c4SBarry Smith 592f4bdf6c4SBarry Smith 593f4bdf6c4SBarry Smith /*MC 594f4bdf6c4SBarry Smith MATHYPRESSTRUCT - MATHYPRESSTRUCT = "hypresstruct" - A matrix type to be used for parallel sparse matrices 595f4bdf6c4SBarry Smith based on the hypre HYPRE_SStructMatrix. 596f4bdf6c4SBarry Smith 597f4bdf6c4SBarry Smith 598f4bdf6c4SBarry Smith Level: intermediate 599f4bdf6c4SBarry Smith 600f4bdf6c4SBarry Smith Notes: Unlike hypre's general semi-struct object consisting of a collection of structured-grid objects and unstructured 601f4bdf6c4SBarry Smith grid objects, we will restrict the semi-struct objects to consist of only structured-grid components. 602f4bdf6c4SBarry Smith 603f4bdf6c4SBarry Smith Unlike the more general support for parts and blocks in hypre this allows only one part, and one block per process and requires the block 604f4bdf6c4SBarry Smith be defined by a DMDA. 605f4bdf6c4SBarry Smith 606c688c046SMatthew G Knepley The matrix needs a DMDA associated with it by either a call to MatSetupDM() or if the matrix is obtained from DMCreateMatrix() 607f4bdf6c4SBarry Smith 608f4bdf6c4SBarry Smith M*/ 609f4bdf6c4SBarry Smith 610f4bdf6c4SBarry Smith #undef __FUNCT__ 611f4bdf6c4SBarry Smith #define __FUNCT__ "MatSetValuesLocal_HYPRESStruct_3d" 612f4bdf6c4SBarry Smith PetscErrorCode MatSetValuesLocal_HYPRESStruct_3d(Mat mat,PetscInt nrow,const PetscInt irow[],PetscInt ncol,const PetscInt icol[],const PetscScalar y[],InsertMode addv) 613f4bdf6c4SBarry Smith { 614f4bdf6c4SBarry Smith PetscErrorCode ierr; 615f4bdf6c4SBarry Smith PetscInt i,j,stencil,index[3]; 616f4bdf6c4SBarry Smith const PetscScalar *values = y; 617f4bdf6c4SBarry Smith Mat_HYPRESStruct *ex = (Mat_HYPRESStruct*) mat->data; 618f4bdf6c4SBarry Smith 6194ddd07fcSJed Brown PetscInt part = 0; /* Petsc sstruct interface only allows 1 part */ 6204ddd07fcSJed Brown PetscInt ordering; 6214ddd07fcSJed Brown PetscInt grid_rank, to_grid_rank; 6224ddd07fcSJed Brown PetscInt var_type, to_var_type; 6234ddd07fcSJed Brown PetscInt to_var_entry = 0; 624f4bdf6c4SBarry Smith 6254ddd07fcSJed Brown PetscInt nvars= ex->nvars; 626f4bdf6c4SBarry Smith PetscInt row,*entries; 627f4bdf6c4SBarry Smith 628f4bdf6c4SBarry Smith PetscFunctionBegin; 629f4bdf6c4SBarry Smith ierr = PetscMalloc(7*nvars*sizeof(PetscInt),&entries);CHKERRQ(ierr); 630f4bdf6c4SBarry Smith ordering= ex->dofs_order; /* ordering= 0 nodal ordering 631f4bdf6c4SBarry Smith 1 variable ordering */ 632f4bdf6c4SBarry Smith /* stencil entries are orderer by variables: var0_stencil0, var0_stencil1, ..., var0_stencil6, var1_stencil0, var1_stencil1, ... */ 633f4bdf6c4SBarry Smith 634f4bdf6c4SBarry Smith if (!ordering) { 635f4bdf6c4SBarry Smith for (i=0; i<nrow; i++) { 636f4bdf6c4SBarry Smith grid_rank= irow[i]/nvars; 637f4bdf6c4SBarry Smith var_type = (irow[i] % nvars); 638f4bdf6c4SBarry Smith 639f4bdf6c4SBarry Smith for (j=0; j<ncol; j++) { 640f4bdf6c4SBarry Smith to_grid_rank= icol[j]/nvars; 641f4bdf6c4SBarry Smith to_var_type = (icol[j] % nvars); 642f4bdf6c4SBarry Smith 643f4bdf6c4SBarry Smith to_var_entry= to_var_entry*7; 644f4bdf6c4SBarry Smith entries[j] = to_var_entry; 645f4bdf6c4SBarry Smith 646f4bdf6c4SBarry Smith stencil = to_grid_rank-grid_rank; 647f4bdf6c4SBarry Smith if (!stencil) { 648f4bdf6c4SBarry Smith entries[j] += 3; 649f4bdf6c4SBarry Smith } else if (stencil == -1) { 650f4bdf6c4SBarry Smith entries[j] += 2; 651f4bdf6c4SBarry Smith } else if (stencil == 1) { 652f4bdf6c4SBarry Smith entries[j] += 4; 653f4bdf6c4SBarry Smith } else if (stencil == -ex->gnx) { 654f4bdf6c4SBarry Smith entries[j] += 1; 655f4bdf6c4SBarry Smith } else if (stencil == ex->gnx) { 656f4bdf6c4SBarry Smith entries[j] += 5; 657f4bdf6c4SBarry Smith } else if (stencil == -ex->gnxgny) { 658f4bdf6c4SBarry Smith entries[j] += 0; 659f4bdf6c4SBarry Smith } else if (stencil == ex->gnxgny) { 660f4bdf6c4SBarry Smith entries[j] += 6; 661f4bdf6c4SBarry Smith } else SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Local row %D local column %D have bad stencil %D",irow[i],icol[j],stencil); 662f4bdf6c4SBarry Smith } 663f4bdf6c4SBarry Smith 664f4bdf6c4SBarry Smith row = ex->gindices[grid_rank] - ex->rstart; 665f4bdf6c4SBarry Smith index[0] = ex->xs + (row % ex->nx); 666f4bdf6c4SBarry Smith index[1] = ex->ys + ((row/ex->nx) % ex->ny); 667f4bdf6c4SBarry Smith index[2] = ex->zs + (row/(ex->nxny)); 668f4bdf6c4SBarry Smith 669f4bdf6c4SBarry Smith if (addv == ADD_VALUES) { 670fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_SStructMatrixAddToValues,(ex->ss_mat,part,index,var_type,ncol,entries,(PetscScalar*)values)); 671f4bdf6c4SBarry Smith } else { 672fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_SStructMatrixSetValues,(ex->ss_mat,part,index,var_type,ncol,entries,(PetscScalar*)values)); 673f4bdf6c4SBarry Smith } 674f4bdf6c4SBarry Smith values += ncol; 675f4bdf6c4SBarry Smith } 676f4bdf6c4SBarry Smith } else { 677f4bdf6c4SBarry Smith for (i=0; i<nrow; i++) { 678f4bdf6c4SBarry Smith var_type = irow[i]/(ex->gnxgnygnz); 679f4bdf6c4SBarry Smith grid_rank= irow[i] - var_type*(ex->gnxgnygnz); 680f4bdf6c4SBarry Smith 681f4bdf6c4SBarry Smith for (j=0; j<ncol; j++) { 682f4bdf6c4SBarry Smith to_var_type = icol[j]/(ex->gnxgnygnz); 683f4bdf6c4SBarry Smith to_grid_rank= icol[j] - to_var_type*(ex->gnxgnygnz); 684f4bdf6c4SBarry Smith 685f4bdf6c4SBarry Smith to_var_entry= to_var_entry*7; 686f4bdf6c4SBarry Smith entries[j] = to_var_entry; 687f4bdf6c4SBarry Smith 688f4bdf6c4SBarry Smith stencil = to_grid_rank-grid_rank; 689f4bdf6c4SBarry Smith if (!stencil) { 690f4bdf6c4SBarry Smith entries[j] += 3; 691f4bdf6c4SBarry Smith } else if (stencil == -1) { 692f4bdf6c4SBarry Smith entries[j] += 2; 693f4bdf6c4SBarry Smith } else if (stencil == 1) { 694f4bdf6c4SBarry Smith entries[j] += 4; 695f4bdf6c4SBarry Smith } else if (stencil == -ex->gnx) { 696f4bdf6c4SBarry Smith entries[j] += 1; 697f4bdf6c4SBarry Smith } else if (stencil == ex->gnx) { 698f4bdf6c4SBarry Smith entries[j] += 5; 699f4bdf6c4SBarry Smith } else if (stencil == -ex->gnxgny) { 700f4bdf6c4SBarry Smith entries[j] += 0; 701f4bdf6c4SBarry Smith } else if (stencil == ex->gnxgny) { 702f4bdf6c4SBarry Smith entries[j] += 6; 703f4bdf6c4SBarry Smith } else SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Local row %D local column %D have bad stencil %D",irow[i],icol[j],stencil); 704f4bdf6c4SBarry Smith } 705f4bdf6c4SBarry Smith 706f4bdf6c4SBarry Smith row = ex->gindices[grid_rank] - ex->rstart; 707f4bdf6c4SBarry Smith index[0] = ex->xs + (row % ex->nx); 708f4bdf6c4SBarry Smith index[1] = ex->ys + ((row/ex->nx) % ex->ny); 709f4bdf6c4SBarry Smith index[2] = ex->zs + (row/(ex->nxny)); 710f4bdf6c4SBarry Smith 711f4bdf6c4SBarry Smith if (addv == ADD_VALUES) { 712fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_SStructMatrixAddToValues,(ex->ss_mat,part,index,var_type,ncol,entries,(PetscScalar*)values)); 713f4bdf6c4SBarry Smith } else { 714fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_SStructMatrixSetValues,(ex->ss_mat,part,index,var_type,ncol,entries,(PetscScalar*)values)); 715f4bdf6c4SBarry Smith } 716f4bdf6c4SBarry Smith values += ncol; 717f4bdf6c4SBarry Smith } 718f4bdf6c4SBarry Smith } 719f4bdf6c4SBarry Smith ierr = PetscFree(entries);CHKERRQ(ierr); 720f4bdf6c4SBarry Smith PetscFunctionReturn(0); 721f4bdf6c4SBarry Smith } 722f4bdf6c4SBarry Smith 723f4bdf6c4SBarry Smith #undef __FUNCT__ 724f4bdf6c4SBarry Smith #define __FUNCT__ "MatZeroRowsLocal_HYPRESStruct_3d" 725f4bdf6c4SBarry Smith PetscErrorCode MatZeroRowsLocal_HYPRESStruct_3d(Mat mat,PetscInt nrow,const PetscInt irow[],PetscScalar d,Vec x,Vec b) 726f4bdf6c4SBarry Smith { 727f4bdf6c4SBarry Smith PetscErrorCode ierr; 728f4bdf6c4SBarry Smith PetscInt i,index[3]; 729f4bdf6c4SBarry Smith PetscScalar **values; 730f4bdf6c4SBarry Smith Mat_HYPRESStruct *ex = (Mat_HYPRESStruct*) mat->data; 731f4bdf6c4SBarry Smith 7324ddd07fcSJed Brown PetscInt part = 0; /* Petsc sstruct interface only allows 1 part */ 7334ddd07fcSJed Brown PetscInt ordering = ex->dofs_order; 7344ddd07fcSJed Brown PetscInt grid_rank; 7354ddd07fcSJed Brown PetscInt var_type; 7364ddd07fcSJed Brown PetscInt nvars= ex->nvars; 737f4bdf6c4SBarry Smith PetscInt row,*entries; 738f4bdf6c4SBarry Smith 739f4bdf6c4SBarry Smith PetscFunctionBegin; 740f4bdf6c4SBarry Smith if (x && b) SETERRQ(((PetscObject)mat)->comm,PETSC_ERR_SUP,"No support"); 741f4bdf6c4SBarry Smith ierr = PetscMalloc(7*nvars*sizeof(PetscInt),&entries);CHKERRQ(ierr); 742f4bdf6c4SBarry Smith 743f4bdf6c4SBarry Smith ierr = PetscMalloc(nvars*sizeof(PetscScalar*),&values);CHKERRQ(ierr); 744f4bdf6c4SBarry Smith ierr = PetscMalloc(7*nvars*nvars*sizeof(PetscScalar),&values[0]);CHKERRQ(ierr); 745f4bdf6c4SBarry Smith for (i=1; i<nvars; i++) { 746f4bdf6c4SBarry Smith values[i] = values[i-1] + nvars*7; 747f4bdf6c4SBarry Smith } 748f4bdf6c4SBarry Smith 749f4bdf6c4SBarry Smith for (i=0; i< nvars; i++) { 750f4bdf6c4SBarry Smith ierr = PetscMemzero(values[i],nvars*7*sizeof(PetscScalar));CHKERRQ(ierr); 751f4bdf6c4SBarry Smith *(values[i]+3) = d; 752f4bdf6c4SBarry Smith } 753f4bdf6c4SBarry Smith 754*8865f1eaSKarl Rupp for (i= 0; i< nvars*7; i++) entries[i] = i; 755f4bdf6c4SBarry Smith 756f4bdf6c4SBarry Smith if (!ordering) { 757f4bdf6c4SBarry Smith for (i=0; i<nrow; i++) { 758f4bdf6c4SBarry Smith grid_rank = irow[i] / nvars; 759f4bdf6c4SBarry Smith var_type =(irow[i] % nvars); 760f4bdf6c4SBarry Smith 761f4bdf6c4SBarry Smith row = ex->gindices[grid_rank] - ex->rstart; 762f4bdf6c4SBarry Smith index[0] = ex->xs + (row % ex->nx); 763f4bdf6c4SBarry Smith index[1] = ex->ys + ((row/ex->nx) % ex->ny); 764f4bdf6c4SBarry Smith index[2] = ex->zs + (row/(ex->nxny)); 765fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_SStructMatrixSetValues,(ex->ss_mat,part,index,var_type,7*nvars,entries,values[var_type])); 766f4bdf6c4SBarry Smith } 767f4bdf6c4SBarry Smith } else { 768f4bdf6c4SBarry Smith for (i=0; i<nrow; i++) { 769f4bdf6c4SBarry Smith var_type = irow[i]/(ex->gnxgnygnz); 770f4bdf6c4SBarry Smith grid_rank= irow[i] - var_type*(ex->gnxgnygnz); 771f4bdf6c4SBarry Smith 772f4bdf6c4SBarry Smith row = ex->gindices[grid_rank] - ex->rstart; 773f4bdf6c4SBarry Smith index[0] = ex->xs + (row % ex->nx); 774f4bdf6c4SBarry Smith index[1] = ex->ys + ((row/ex->nx) % ex->ny); 775f4bdf6c4SBarry Smith index[2] = ex->zs + (row/(ex->nxny)); 776fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_SStructMatrixSetValues,(ex->ss_mat,part,index,var_type,7*nvars,entries,values[var_type])); 777f4bdf6c4SBarry Smith } 778f4bdf6c4SBarry Smith } 779fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_SStructMatrixAssemble,(ex->ss_mat)); 780f4bdf6c4SBarry Smith ierr = PetscFree(values[0]);CHKERRQ(ierr); 781f4bdf6c4SBarry Smith ierr = PetscFree(values);CHKERRQ(ierr); 782f4bdf6c4SBarry Smith ierr = PetscFree(entries);CHKERRQ(ierr); 783f4bdf6c4SBarry Smith PetscFunctionReturn(0); 784f4bdf6c4SBarry Smith } 785f4bdf6c4SBarry Smith 786f4bdf6c4SBarry Smith #undef __FUNCT__ 787f4bdf6c4SBarry Smith #define __FUNCT__ "MatZeroEntries_HYPRESStruct_3d" 788f4bdf6c4SBarry Smith PetscErrorCode MatZeroEntries_HYPRESStruct_3d(Mat mat) 789f4bdf6c4SBarry Smith { 790f4bdf6c4SBarry Smith PetscErrorCode ierr; 791f4bdf6c4SBarry Smith Mat_HYPRESStruct *ex = (Mat_HYPRESStruct*) mat->data; 7924ddd07fcSJed Brown PetscInt nvars= ex->nvars; 7934ddd07fcSJed Brown PetscInt size; 7944ddd07fcSJed Brown PetscInt part= 0; /* only one part */ 795f4bdf6c4SBarry Smith 796f4bdf6c4SBarry Smith PetscFunctionBegin; 797f4bdf6c4SBarry Smith size = ((ex->hbox.imax[0])-(ex->hbox.imin[0])+1)*((ex->hbox.imax[1])-(ex->hbox.imin[1])+1)*((ex->hbox.imax[2])-(ex->hbox.imin[2])+1); 798f4bdf6c4SBarry Smith { 7994ddd07fcSJed Brown PetscInt i,*entries,iupper[3],ilower[3]; 800f4bdf6c4SBarry Smith PetscScalar *values; 8014ddd07fcSJed Brown 802f4bdf6c4SBarry Smith 803f4bdf6c4SBarry Smith for (i= 0; i< 3; i++) { 804f4bdf6c4SBarry Smith ilower[i]= ex->hbox.imin[i]; 805f4bdf6c4SBarry Smith iupper[i]= ex->hbox.imax[i]; 806f4bdf6c4SBarry Smith } 807f4bdf6c4SBarry Smith 808f4bdf6c4SBarry Smith ierr = PetscMalloc2(nvars*7,PetscInt,&entries,nvars*7*size,PetscScalar,&values);CHKERRQ(ierr); 809*8865f1eaSKarl Rupp for (i= 0; i< nvars*7; i++) entries[i]= i; 810f4bdf6c4SBarry Smith ierr = PetscMemzero(values,nvars*7*size*sizeof(PetscScalar));CHKERRQ(ierr); 811f4bdf6c4SBarry Smith 812f4bdf6c4SBarry Smith for (i= 0; i< nvars; i++) { 813fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_SStructMatrixSetBoxValues,(ex->ss_mat,part,ilower,iupper,i,nvars*7,entries,values)); 814f4bdf6c4SBarry Smith } 815f4bdf6c4SBarry Smith ierr = PetscFree2(entries,values);CHKERRQ(ierr); 816f4bdf6c4SBarry Smith } 817fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_SStructMatrixAssemble,(ex->ss_mat)); 818f4bdf6c4SBarry Smith PetscFunctionReturn(0); 819f4bdf6c4SBarry Smith } 820f4bdf6c4SBarry Smith 821f4bdf6c4SBarry Smith 822f4bdf6c4SBarry Smith #undef __FUNCT__ 823c688c046SMatthew G Knepley #define __FUNCT__ "MatSetupDM_HYPRESStruct" 824c688c046SMatthew G Knepley PetscErrorCode MatSetupDM_HYPRESStruct(Mat mat,DM da) 825f4bdf6c4SBarry Smith { 826f4bdf6c4SBarry Smith PetscErrorCode ierr; 827f4bdf6c4SBarry Smith Mat_HYPRESStruct *ex = (Mat_HYPRESStruct*) mat->data; 828f4bdf6c4SBarry Smith PetscInt dim,dof,sw[3],nx,ny,nz; 8294ddd07fcSJed Brown PetscInt ilower[3],iupper[3],ssize,i; 8307ab44359SJed Brown DMDABoundaryType px,py,pz; 831f4bdf6c4SBarry Smith DMDAStencilType st; 8324ddd07fcSJed Brown PetscInt nparts= 1; /* assuming only one part */ 8334ddd07fcSJed Brown PetscInt part = 0; 834f4bdf6c4SBarry Smith 835f4bdf6c4SBarry Smith PetscFunctionBegin; 836f4bdf6c4SBarry Smith ex->da = da; 837f4bdf6c4SBarry Smith ierr = PetscObjectReference((PetscObject)da);CHKERRQ(ierr); 838f4bdf6c4SBarry Smith 8397ab44359SJed Brown ierr = DMDAGetInfo(ex->da,&dim,0,0,0,0,0,0,&dof,&sw[0],&px,&py,&pz,&st);CHKERRQ(ierr); 840f4bdf6c4SBarry Smith ierr = DMDAGetCorners(ex->da,&ilower[0],&ilower[1],&ilower[2],&iupper[0],&iupper[1],&iupper[2]);CHKERRQ(ierr); 841f4bdf6c4SBarry Smith iupper[0] += ilower[0] - 1; 842f4bdf6c4SBarry Smith iupper[1] += ilower[1] - 1; 843f4bdf6c4SBarry Smith iupper[2] += ilower[2] - 1; 844f4bdf6c4SBarry Smith /* the hypre_Box is used to zero out the matrix entries in MatZeroValues() */ 845f4bdf6c4SBarry Smith ex->hbox.imin[0] = ilower[0]; 846f4bdf6c4SBarry Smith ex->hbox.imin[1] = ilower[1]; 847f4bdf6c4SBarry Smith ex->hbox.imin[2] = ilower[2]; 848f4bdf6c4SBarry Smith ex->hbox.imax[0] = iupper[0]; 849f4bdf6c4SBarry Smith ex->hbox.imax[1] = iupper[1]; 850f4bdf6c4SBarry Smith ex->hbox.imax[2] = iupper[2]; 851f4bdf6c4SBarry Smith 852f4bdf6c4SBarry Smith ex->dofs_order = 0; 853f4bdf6c4SBarry Smith 854f4bdf6c4SBarry Smith /* assuming that the same number of dofs on each gridpoint. Also assume all cell-centred based */ 855f4bdf6c4SBarry Smith ex->nvars= dof; 856f4bdf6c4SBarry Smith 857f4bdf6c4SBarry Smith /* create the hypre grid object and set its information */ 8587ab44359SJed Brown if (px || py || pz) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"Ask us to add periodic support by calling HYPRE_SStructGridSetPeriodic()"); 859fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_SStructGridCreate,(ex->hcomm,dim,nparts,&ex->ss_grid)); 860f4bdf6c4SBarry Smith 861fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_SStructGridSetExtents,(ex->ss_grid,part,ex->hbox.imin,ex->hbox.imax)); 862f4bdf6c4SBarry Smith 863f4bdf6c4SBarry Smith { 864f4bdf6c4SBarry Smith HYPRE_SStructVariable *vartypes; 865f4bdf6c4SBarry Smith ierr = PetscMalloc(ex->nvars*sizeof(HYPRE_SStructVariable),&vartypes);CHKERRQ(ierr); 866*8865f1eaSKarl Rupp for (i= 0; i< ex->nvars; i++) vartypes[i]= HYPRE_SSTRUCT_VARIABLE_CELL; 867fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_SStructGridSetVariables,(ex->ss_grid, part, ex->nvars,vartypes)); 868f4bdf6c4SBarry Smith ierr = PetscFree(vartypes);CHKERRQ(ierr); 869f4bdf6c4SBarry Smith } 870fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_SStructGridAssemble,(ex->ss_grid)); 871f4bdf6c4SBarry Smith 872f4bdf6c4SBarry Smith sw[1] = sw[0]; 873f4bdf6c4SBarry Smith sw[2] = sw[1]; 874fd3f9acdSBarry Smith /* PetscStackCallStandard(HYPRE_SStructGridSetNumGhost,(ex->ss_grid,sw)); */ 875f4bdf6c4SBarry Smith 876f4bdf6c4SBarry Smith /* create the hypre stencil object and set its information */ 877f4bdf6c4SBarry Smith if (sw[0] > 1) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"Ask us to add support for wider stencils"); 878f4bdf6c4SBarry Smith if (st == DMDA_STENCIL_BOX) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"Ask us to add support for box stencils"); 879f4bdf6c4SBarry Smith 880f4bdf6c4SBarry Smith if (dim == 1) { 8814ddd07fcSJed Brown PetscInt offsets[3][1] = {{-1},{0},{1}}; 8824ddd07fcSJed Brown PetscInt j, cnt; 883f4bdf6c4SBarry Smith 884f4bdf6c4SBarry Smith ssize = 3*(ex->nvars); 885fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_SStructStencilCreate,(dim,ssize,&ex->ss_stencil)); 886f4bdf6c4SBarry Smith cnt= 0; 887f4bdf6c4SBarry Smith for (i = 0; i < (ex->nvars); i++) { 888f4bdf6c4SBarry Smith for (j = 0; j < 3; j++) { 889fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_SStructStencilSetEntry,(ex->ss_stencil, cnt, offsets[j], i)); 890f4bdf6c4SBarry Smith cnt++; 891f4bdf6c4SBarry Smith } 892f4bdf6c4SBarry Smith } 893f4bdf6c4SBarry Smith } else if (dim == 2) { 8944ddd07fcSJed Brown PetscInt offsets[5][2] = {{0,-1},{-1,0},{0,0},{1,0},{0,1}}; 8954ddd07fcSJed Brown PetscInt j, cnt; 896f4bdf6c4SBarry Smith 897f4bdf6c4SBarry Smith ssize = 5*(ex->nvars); 898fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_SStructStencilCreate,(dim,ssize,&ex->ss_stencil)); 899f4bdf6c4SBarry Smith cnt= 0; 900f4bdf6c4SBarry Smith for (i= 0; i< (ex->nvars); i++) { 901f4bdf6c4SBarry Smith for (j= 0; j< 5; j++) { 902fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_SStructStencilSetEntry,(ex->ss_stencil, cnt, offsets[j], i)); 903f4bdf6c4SBarry Smith cnt++; 904f4bdf6c4SBarry Smith } 905f4bdf6c4SBarry Smith } 906f4bdf6c4SBarry Smith } else if (dim == 3) { 9074ddd07fcSJed Brown PetscInt offsets[7][3] = {{0,0,-1},{0,-1,0},{-1,0,0},{0,0,0},{1,0,0},{0,1,0},{0,0,1}}; 9084ddd07fcSJed Brown PetscInt j, cnt; 909f4bdf6c4SBarry Smith 910f4bdf6c4SBarry Smith ssize = 7*(ex->nvars); 911fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_SStructStencilCreate,(dim,ssize,&ex->ss_stencil)); 912f4bdf6c4SBarry Smith cnt= 0; 913f4bdf6c4SBarry Smith for (i= 0; i< (ex->nvars); i++) { 914f4bdf6c4SBarry Smith for (j= 0; j< 7; j++) { 915fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_SStructStencilSetEntry,(ex->ss_stencil, cnt, offsets[j], i)); 916f4bdf6c4SBarry Smith cnt++; 917f4bdf6c4SBarry Smith } 918f4bdf6c4SBarry Smith } 919f4bdf6c4SBarry Smith } 920f4bdf6c4SBarry Smith 921f4bdf6c4SBarry Smith /* create the HYPRE graph */ 922fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_SStructGraphCreate,(ex->hcomm, ex->ss_grid, &(ex->ss_graph))); 923f4bdf6c4SBarry Smith 924f4bdf6c4SBarry Smith /* set the stencil graph. Note that each variable has the same graph. This means that each 925f4bdf6c4SBarry Smith variable couples to all the other variable and with the same stencil pattern. */ 926f4bdf6c4SBarry Smith for (i= 0; i< (ex->nvars); i++) { 927fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_SStructGraphSetStencil,(ex->ss_graph,part,i,ex->ss_stencil)); 928f4bdf6c4SBarry Smith } 929fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_SStructGraphAssemble,(ex->ss_graph)); 930f4bdf6c4SBarry Smith 931f4bdf6c4SBarry Smith /* create the HYPRE sstruct vectors for rhs and solution */ 932fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_SStructVectorCreate,(ex->hcomm,ex->ss_grid,&ex->ss_b)); 933fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_SStructVectorCreate,(ex->hcomm,ex->ss_grid,&ex->ss_x)); 934fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_SStructVectorInitialize,(ex->ss_b)); 935fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_SStructVectorInitialize,(ex->ss_x)); 936fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_SStructVectorAssemble,(ex->ss_b)); 937fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_SStructVectorAssemble,(ex->ss_x)); 938f4bdf6c4SBarry Smith 939f4bdf6c4SBarry Smith /* create the hypre matrix object and set its information */ 940fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_SStructMatrixCreate,(ex->hcomm,ex->ss_graph,&ex->ss_mat)); 941fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_SStructGridDestroy,(ex->ss_grid)); 942fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_SStructStencilDestroy,(ex->ss_stencil)); 943f4bdf6c4SBarry Smith if (ex->needsinitialization) { 944fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_SStructMatrixInitialize,(ex->ss_mat)); 945f4bdf6c4SBarry Smith ex->needsinitialization = PETSC_FALSE; 946f4bdf6c4SBarry Smith } 947f4bdf6c4SBarry Smith 948f4bdf6c4SBarry Smith /* set the global and local sizes of the matrix */ 949f4bdf6c4SBarry Smith ierr = DMDAGetCorners(da,0,0,0,&nx,&ny,&nz);CHKERRQ(ierr); 950f4bdf6c4SBarry Smith ierr = MatSetSizes(mat,dof*nx*ny*nz,dof*nx*ny*nz,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr); 951f4bdf6c4SBarry Smith ierr = PetscLayoutSetBlockSize(mat->rmap,1);CHKERRQ(ierr); 952f4bdf6c4SBarry Smith ierr = PetscLayoutSetBlockSize(mat->cmap,1);CHKERRQ(ierr); 953f4bdf6c4SBarry Smith ierr = PetscLayoutSetUp(mat->rmap);CHKERRQ(ierr); 954f4bdf6c4SBarry Smith ierr = PetscLayoutSetUp(mat->cmap);CHKERRQ(ierr); 955f4bdf6c4SBarry Smith 956f4bdf6c4SBarry Smith if (dim == 3) { 957f4bdf6c4SBarry Smith mat->ops->setvalueslocal = MatSetValuesLocal_HYPRESStruct_3d; 958f4bdf6c4SBarry Smith mat->ops->zerorowslocal = MatZeroRowsLocal_HYPRESStruct_3d; 959f4bdf6c4SBarry Smith mat->ops->zeroentries = MatZeroEntries_HYPRESStruct_3d; 960*8865f1eaSKarl Rupp 961f4bdf6c4SBarry Smith ierr = MatZeroEntries_HYPRESStruct_3d(mat);CHKERRQ(ierr); 962f4bdf6c4SBarry Smith } else SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"Only support for 3d DMDA currently"); 963f4bdf6c4SBarry Smith 964f4bdf6c4SBarry Smith /* get values that will be used repeatedly in MatSetValuesLocal() and MatZeroRowsLocal() repeatedly */ 965f4bdf6c4SBarry Smith ierr = MatGetOwnershipRange(mat,&ex->rstart,PETSC_NULL);CHKERRQ(ierr); 966f4bdf6c4SBarry Smith ierr = DMDAGetGlobalIndices(ex->da,PETSC_NULL,&ex->gindices);CHKERRQ(ierr); 967f4bdf6c4SBarry Smith ierr = DMDAGetGhostCorners(ex->da,0,0,0,&ex->gnx,&ex->gnxgny,&ex->gnxgnygnz);CHKERRQ(ierr); 968*8865f1eaSKarl Rupp 969f4bdf6c4SBarry Smith ex->gnxgny *= ex->gnx; 970f4bdf6c4SBarry Smith ex->gnxgnygnz *= ex->gnxgny; 971*8865f1eaSKarl Rupp 972f4bdf6c4SBarry Smith ierr = DMDAGetCorners(ex->da,&ex->xs,&ex->ys,&ex->zs,&ex->nx,&ex->ny,&ex->nz);CHKERRQ(ierr); 973*8865f1eaSKarl Rupp 974f4bdf6c4SBarry Smith ex->nxny = ex->nx*ex->ny; 975f4bdf6c4SBarry Smith ex->nxnynz = ex->nz*ex->nxny; 976f4bdf6c4SBarry Smith PetscFunctionReturn(0); 977f4bdf6c4SBarry Smith } 978f4bdf6c4SBarry Smith 979f4bdf6c4SBarry Smith #undef __FUNCT__ 980f4bdf6c4SBarry Smith #define __FUNCT__ "MatMult_HYPRESStruct" 981f4bdf6c4SBarry Smith PetscErrorCode MatMult_HYPRESStruct(Mat A,Vec x,Vec y) 982f4bdf6c4SBarry Smith { 983f4bdf6c4SBarry Smith PetscErrorCode ierr; 984f4bdf6c4SBarry Smith PetscScalar *xx,*yy; 9854ddd07fcSJed Brown PetscInt ilower[3],iupper[3]; 986f4bdf6c4SBarry Smith Mat_HYPRESStruct *mx = (Mat_HYPRESStruct*)(A->data); 9874ddd07fcSJed Brown PetscInt ordering= mx->dofs_order; 9884ddd07fcSJed Brown PetscInt nvars = mx->nvars; 9894ddd07fcSJed Brown PetscInt part = 0; 9904ddd07fcSJed Brown PetscInt size; 9914ddd07fcSJed Brown PetscInt i; 992f4bdf6c4SBarry Smith 993f4bdf6c4SBarry Smith PetscFunctionBegin; 994f4bdf6c4SBarry Smith ierr = DMDAGetCorners(mx->da,&ilower[0],&ilower[1],&ilower[2],&iupper[0],&iupper[1],&iupper[2]);CHKERRQ(ierr); 995f4bdf6c4SBarry Smith iupper[0] += ilower[0] - 1; 996f4bdf6c4SBarry Smith iupper[1] += ilower[1] - 1; 997f4bdf6c4SBarry Smith iupper[2] += ilower[2] - 1; 998f4bdf6c4SBarry Smith 999f4bdf6c4SBarry Smith size= 1; 1000*8865f1eaSKarl Rupp for (i = 0; i < 3; i++) size *= (iupper[i]-ilower[i]+1); 1001f4bdf6c4SBarry Smith 1002f4bdf6c4SBarry Smith /* copy x values over to hypre for variable ordering */ 1003f4bdf6c4SBarry Smith if (ordering) { 1004fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_SStructVectorSetConstantValues,(mx->ss_b,0.0)); 1005f4bdf6c4SBarry Smith ierr = VecGetArray(x,&xx);CHKERRQ(ierr); 1006f4bdf6c4SBarry Smith for (i= 0; i< nvars; i++) { 1007fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_SStructVectorSetBoxValues,(mx->ss_b,part,ilower,iupper,i,xx+(size*i))); 1008f4bdf6c4SBarry Smith } 1009f4bdf6c4SBarry Smith ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr); 1010fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_SStructVectorAssemble,(mx->ss_b)); 1011fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_SStructMatrixMatvec,(1.0,mx->ss_mat,mx->ss_b,0.0,mx->ss_x)); 1012f4bdf6c4SBarry Smith 1013f4bdf6c4SBarry Smith /* copy solution values back to PETSc */ 1014f4bdf6c4SBarry Smith ierr = VecGetArray(y,&yy);CHKERRQ(ierr); 1015f4bdf6c4SBarry Smith for (i= 0; i< nvars; i++) { 1016fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_SStructVectorGetBoxValues,(mx->ss_x,part,ilower,iupper,i,yy+(size*i))); 1017f4bdf6c4SBarry Smith } 1018f4bdf6c4SBarry Smith ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr); 1019f4bdf6c4SBarry Smith } else { /* nodal ordering must be mapped to variable ordering for sys_pfmg */ 1020f4bdf6c4SBarry Smith PetscScalar *z; 10214ddd07fcSJed Brown PetscInt j, k; 1022f4bdf6c4SBarry Smith 1023f4bdf6c4SBarry Smith ierr = PetscMalloc(nvars*size*sizeof(PetscScalar),&z);CHKERRQ(ierr); 1024fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_SStructVectorSetConstantValues,(mx->ss_b,0.0)); 1025f4bdf6c4SBarry Smith ierr = VecGetArray(x,&xx);CHKERRQ(ierr); 1026f4bdf6c4SBarry Smith 1027f4bdf6c4SBarry Smith /* transform nodal to hypre's variable ordering for sys_pfmg */ 1028f4bdf6c4SBarry Smith for (i= 0; i< size; i++) { 1029f4bdf6c4SBarry Smith k= i*nvars; 1030*8865f1eaSKarl Rupp for (j= 0; j< nvars; j++) z[j*size+i]= xx[k+j]; 1031f4bdf6c4SBarry Smith } 1032f4bdf6c4SBarry Smith for (i= 0; i< nvars; i++) { 1033fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_SStructVectorSetBoxValues,(mx->ss_b,part,ilower,iupper,i,z+(size*i))); 1034f4bdf6c4SBarry Smith } 1035f4bdf6c4SBarry Smith ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr); 1036fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_SStructVectorAssemble,(mx->ss_b)); 1037fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_SStructMatrixMatvec,(1.0,mx->ss_mat,mx->ss_b,0.0,mx->ss_x)); 1038f4bdf6c4SBarry Smith 1039f4bdf6c4SBarry Smith /* copy solution values back to PETSc */ 1040f4bdf6c4SBarry Smith ierr = VecGetArray(y,&yy);CHKERRQ(ierr); 1041f4bdf6c4SBarry Smith for (i= 0; i< nvars; i++) { 1042fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_SStructVectorGetBoxValues,(mx->ss_x,part,ilower,iupper,i,z+(size*i))); 1043f4bdf6c4SBarry Smith } 1044f4bdf6c4SBarry Smith /* transform hypre's variable ordering for sys_pfmg to nodal ordering */ 1045f4bdf6c4SBarry Smith for (i= 0; i< size; i++) { 1046f4bdf6c4SBarry Smith k= i*nvars; 1047*8865f1eaSKarl Rupp for (j= 0; j< nvars; j++) yy[k+j]= z[j*size+i]; 1048f4bdf6c4SBarry Smith } 1049f4bdf6c4SBarry Smith ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr); 1050f4bdf6c4SBarry Smith ierr = PetscFree(z);CHKERRQ(ierr); 1051f4bdf6c4SBarry Smith } 1052f4bdf6c4SBarry Smith PetscFunctionReturn(0); 1053f4bdf6c4SBarry Smith } 1054f4bdf6c4SBarry Smith 1055f4bdf6c4SBarry Smith #undef __FUNCT__ 1056f4bdf6c4SBarry Smith #define __FUNCT__ "MatAssemblyEnd_HYPRESStruct" 1057f4bdf6c4SBarry Smith PetscErrorCode MatAssemblyEnd_HYPRESStruct(Mat mat,MatAssemblyType mode) 1058f4bdf6c4SBarry Smith { 1059f4bdf6c4SBarry Smith Mat_HYPRESStruct *ex = (Mat_HYPRESStruct*) mat->data; 1060f4bdf6c4SBarry Smith PetscErrorCode ierr; 1061f4bdf6c4SBarry Smith 1062f4bdf6c4SBarry Smith PetscFunctionBegin; 1063fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_SStructMatrixAssemble,(ex->ss_mat)); 1064f4bdf6c4SBarry Smith PetscFunctionReturn(0); 1065f4bdf6c4SBarry Smith } 1066f4bdf6c4SBarry Smith 1067f4bdf6c4SBarry Smith #undef __FUNCT__ 1068f4bdf6c4SBarry Smith #define __FUNCT__ "MatZeroEntries_HYPRESStruct" 1069f4bdf6c4SBarry Smith PetscErrorCode MatZeroEntries_HYPRESStruct(Mat mat) 1070f4bdf6c4SBarry Smith { 1071f4bdf6c4SBarry Smith PetscFunctionBegin; 1072f4bdf6c4SBarry Smith /* before the DMDA is set to the matrix the zero doesn't need to do anything */ 1073f4bdf6c4SBarry Smith PetscFunctionReturn(0); 1074f4bdf6c4SBarry Smith } 1075f4bdf6c4SBarry Smith 1076f4bdf6c4SBarry Smith 1077f4bdf6c4SBarry Smith #undef __FUNCT__ 1078f4bdf6c4SBarry Smith #define __FUNCT__ "MatDestroy_HYPRESStruct" 1079f4bdf6c4SBarry Smith PetscErrorCode MatDestroy_HYPRESStruct(Mat mat) 1080f4bdf6c4SBarry Smith { 1081f4bdf6c4SBarry Smith Mat_HYPRESStruct *ex = (Mat_HYPRESStruct*) mat->data; 1082f4bdf6c4SBarry Smith PetscErrorCode ierr; 1083f4bdf6c4SBarry Smith 1084f4bdf6c4SBarry Smith PetscFunctionBegin; 1085fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_SStructGraphDestroy,(ex->ss_graph)); 1086fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_SStructMatrixDestroy,(ex->ss_mat)); 1087fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_SStructVectorDestroy,(ex->ss_x)); 1088fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_SStructVectorDestroy,(ex->ss_b)); 1089f4bdf6c4SBarry Smith PetscFunctionReturn(0); 1090f4bdf6c4SBarry Smith } 1091f4bdf6c4SBarry Smith 1092f4bdf6c4SBarry Smith EXTERN_C_BEGIN 1093f4bdf6c4SBarry Smith #undef __FUNCT__ 1094f4bdf6c4SBarry Smith #define __FUNCT__ "MatCreate_HYPRESStruct" 1095f4bdf6c4SBarry Smith PetscErrorCode MatCreate_HYPRESStruct(Mat B) 1096f4bdf6c4SBarry Smith { 1097f4bdf6c4SBarry Smith Mat_HYPRESStruct *ex; 1098f4bdf6c4SBarry Smith PetscErrorCode ierr; 1099f4bdf6c4SBarry Smith 1100f4bdf6c4SBarry Smith PetscFunctionBegin; 1101f4bdf6c4SBarry Smith ierr = PetscNewLog(B,Mat_HYPRESStruct,&ex);CHKERRQ(ierr); 1102f4bdf6c4SBarry Smith B->data = (void*)ex; 1103f4bdf6c4SBarry Smith B->rmap->bs = 1; 1104f4bdf6c4SBarry Smith B->assembled = PETSC_FALSE; 1105f4bdf6c4SBarry Smith 1106f4bdf6c4SBarry Smith B->insertmode = NOT_SET_VALUES; 1107f4bdf6c4SBarry Smith 1108f4bdf6c4SBarry Smith B->ops->assemblyend = MatAssemblyEnd_HYPRESStruct; 1109f4bdf6c4SBarry Smith B->ops->mult = MatMult_HYPRESStruct; 1110f4bdf6c4SBarry Smith B->ops->zeroentries = MatZeroEntries_HYPRESStruct; 1111f4bdf6c4SBarry Smith B->ops->destroy = MatDestroy_HYPRESStruct; 1112f4bdf6c4SBarry Smith 1113f4bdf6c4SBarry Smith ex->needsinitialization = PETSC_TRUE; 1114f4bdf6c4SBarry Smith 1115f4bdf6c4SBarry Smith ierr = MPI_Comm_dup(((PetscObject)B)->comm,&(ex->hcomm));CHKERRQ(ierr); 1116c688c046SMatthew G Knepley ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSetupDM_C","MatSetupDM_HYPRESStruct",MatSetupDM_HYPRESStruct);CHKERRQ(ierr); 1117f4bdf6c4SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)B,MATHYPRESSTRUCT);CHKERRQ(ierr); 1118f4bdf6c4SBarry Smith PetscFunctionReturn(0); 1119f4bdf6c4SBarry Smith } 1120f4bdf6c4SBarry Smith EXTERN_C_END 1121f4bdf6c4SBarry Smith 1122f4bdf6c4SBarry Smith 1123f4bdf6c4SBarry Smith 1124