1f4bdf6c4SBarry Smith 2f4bdf6c4SBarry Smith /* 3f4bdf6c4SBarry Smith Creates hypre ijmatrix from PETSc matrix 4f4bdf6c4SBarry Smith */ 5c6db04a5SJed Brown #include <petscsys.h> 6c6db04a5SJed Brown #include <private/matimpl.h> /*I "petscmat.h" I*/ 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; 15f4bdf6c4SBarry Smith PetscInt i; 16f4bdf6c4SBarry Smith PetscInt n_d,*ia_d,n_o,*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 } 48f4bdf6c4SBarry Smith PetscStackCallHypre(0,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; 60f4bdf6c4SBarry Smith int 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); 66*4994cf47SJed 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; 71f4bdf6c4SBarry Smith PetscStackCallHypre(0,HYPRE_IJMatrixCreate,(((PetscObject)A)->comm,rstart,rend-1,cstart,cend-1,ij)); 72f4bdf6c4SBarry Smith PetscStackCallHypre(0,HYPRE_IJMatrixSetObjectType,(*ij,HYPRE_PARCSR)); 73f4bdf6c4SBarry Smith { 74f4bdf6c4SBarry Smith PetscBool same; 75f4bdf6c4SBarry Smith Mat A_d,A_o; 76f4bdf6c4SBarry Smith PetscInt *colmap; 77f4bdf6c4SBarry Smith ierr = PetscTypeCompare((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 } 83f4bdf6c4SBarry Smith ierr = PetscTypeCompare((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 } 89f4bdf6c4SBarry Smith ierr = PetscTypeCompare((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 } 94f4bdf6c4SBarry Smith ierr = PetscTypeCompare((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; 120f4bdf6c4SBarry Smith ierr = PetscTypeCompare((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 } 125f4bdf6c4SBarry Smith ierr = PetscTypeCompare((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); 132f4bdf6c4SBarry Smith PetscStackCallHypre(0,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); 136f4bdf6c4SBarry Smith PetscStackCallHypre(0,HYPRE_IJMatrixSetValues,(ij,1,&ncols,&i,cols,values)); 137f4bdf6c4SBarry Smith ierr = MatRestoreRow(A,i,&ncols,&cols,&values);CHKERRQ(ierr); 138f4bdf6c4SBarry Smith } 139f4bdf6c4SBarry Smith PetscStackCallHypre(0,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); 170f4bdf6c4SBarry Smith PetscStackCallHypre(0,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; 183f4bdf6c4SBarry Smith PetscStackCallHypre(0,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); 211f4bdf6c4SBarry Smith PetscStackCallHypre(0,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; 224f4bdf6c4SBarry Smith for (i=0; i<pdiag->nz; i++) { 225f4bdf6c4SBarry Smith jj[i] = cstart + pjj[i]; 226f4bdf6c4SBarry Smith } 227f4bdf6c4SBarry Smith ierr = PetscMemcpy(hdiag->data,pdiag->a,pdiag->nz*sizeof(PetscScalar)); 228f4bdf6c4SBarry Smith ierr = PetscMemcpy(hoffd->i,poffd->i,(pA->A->rmap->n + 1)*sizeof(PetscInt)); 229f4bdf6c4SBarry Smith /* need to move the offd column indices (hoffd->j) back to global numbering since hypre is expecting this 230f4bdf6c4SBarry Smith If we hacked a hypre a bit more we might be able to avoid this step */ 231f4bdf6c4SBarry Smith jj = hoffd->j; 232f4bdf6c4SBarry Smith pjj = poffd->j; 233f4bdf6c4SBarry Smith for (i=0; i<poffd->nz; i++) { 234f4bdf6c4SBarry Smith jj[i] = garray[pjj[i]]; 235f4bdf6c4SBarry Smith } 236f4bdf6c4SBarry Smith ierr = PetscMemcpy(hoffd->data,poffd->a,poffd->nz*sizeof(PetscScalar)); 237f4bdf6c4SBarry Smith 238f4bdf6c4SBarry Smith hypre_AuxParCSRMatrixNeedAux(aux_matrix) = 0; 239f4bdf6c4SBarry Smith PetscStackCallHypre(0,HYPRE_IJMatrixAssemble,(ij)); 240f4bdf6c4SBarry Smith ierr = PetscLogEventEnd(MAT_Convert,A,0,0,0);CHKERRQ(ierr); 241f4bdf6c4SBarry Smith PetscFunctionReturn(0); 242f4bdf6c4SBarry Smith } 243f4bdf6c4SBarry Smith 244f4bdf6c4SBarry Smith /* 245f4bdf6c4SBarry Smith Does NOT copy the data over, instead uses DIRECTLY the pointers from the PETSc MPIAIJ format 246f4bdf6c4SBarry Smith 247f4bdf6c4SBarry Smith This is UNFINISHED and does NOT work! The problem is that hypre puts the diagonal entry first 248f4bdf6c4SBarry Smith which will corrupt the PETSc data structure if we did this. Need a work around to this problem. 249f4bdf6c4SBarry Smith */ 250c6db04a5SJed Brown #include <_hypre_IJ_mv.h> 251c6db04a5SJed Brown #include <HYPRE_IJ_mv.h> 252f4bdf6c4SBarry Smith 253f4bdf6c4SBarry Smith #undef __FUNCT__ 254f4bdf6c4SBarry Smith #define __FUNCT__ "MatHYPRE_IJMatrixLink" 255f4bdf6c4SBarry Smith PetscErrorCode MatHYPRE_IJMatrixLink(Mat A,HYPRE_IJMatrix *ij) 256f4bdf6c4SBarry Smith { 257f4bdf6c4SBarry Smith PetscErrorCode ierr; 258f4bdf6c4SBarry Smith int rstart,rend,cstart,cend; 259f4bdf6c4SBarry Smith PetscBool flg; 260f4bdf6c4SBarry Smith hypre_AuxParCSRMatrix *aux_matrix; 261f4bdf6c4SBarry Smith 262f4bdf6c4SBarry Smith PetscFunctionBegin; 263f4bdf6c4SBarry Smith PetscValidHeaderSpecific(A,MAT_CLASSID,1); 264f4bdf6c4SBarry Smith PetscValidType(A,1); 265f4bdf6c4SBarry Smith PetscValidPointer(ij,2); 266f4bdf6c4SBarry Smith ierr = PetscTypeCompare((PetscObject)A,MATMPIAIJ,&flg);CHKERRQ(ierr); 267f4bdf6c4SBarry Smith if (!flg) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"Can only use with PETSc MPIAIJ matrices"); 268*4994cf47SJed Brown ierr = MatSetUp(A);CHKERRQ(ierr); 269f4bdf6c4SBarry Smith 270f4bdf6c4SBarry Smith ierr = PetscLogEventBegin(MAT_Convert,A,0,0,0);CHKERRQ(ierr); 271f4bdf6c4SBarry Smith rstart = A->rmap->rstart; 272f4bdf6c4SBarry Smith rend = A->rmap->rend; 273f4bdf6c4SBarry Smith cstart = A->cmap->rstart; 274f4bdf6c4SBarry Smith cend = A->cmap->rend; 275f4bdf6c4SBarry Smith PetscStackCallHypre(0,HYPRE_IJMatrixCreate,(((PetscObject)A)->comm,rstart,rend-1,cstart,cend-1,ij)); 276f4bdf6c4SBarry Smith PetscStackCallHypre(0,HYPRE_IJMatrixSetObjectType,(*ij,HYPRE_PARCSR)); 277f4bdf6c4SBarry Smith 278f4bdf6c4SBarry Smith PetscStackCallHypre(0,HYPRE_IJMatrixInitialize,(*ij)); 279f4bdf6c4SBarry Smith aux_matrix = (hypre_AuxParCSRMatrix*)hypre_IJMatrixTranslator(*ij); 280f4bdf6c4SBarry Smith 281f4bdf6c4SBarry Smith hypre_AuxParCSRMatrixNeedAux(aux_matrix) = 0; 282f4bdf6c4SBarry Smith 283f4bdf6c4SBarry Smith /* this is the Hack part where we monkey directly with the hypre datastructures */ 284f4bdf6c4SBarry Smith 285f4bdf6c4SBarry Smith PetscStackCallHypre(0,HYPRE_IJMatrixAssemble,(*ij)); 286f4bdf6c4SBarry Smith ierr = PetscLogEventEnd(MAT_Convert,A,0,0,0);CHKERRQ(ierr); 287f4bdf6c4SBarry Smith PetscFunctionReturn(0); 288f4bdf6c4SBarry Smith } 289f4bdf6c4SBarry Smith 290f4bdf6c4SBarry Smith /* -----------------------------------------------------------------------------------------------------------------*/ 291f4bdf6c4SBarry Smith 292f4bdf6c4SBarry Smith /*MC 293f4bdf6c4SBarry Smith MATHYPRESTRUCT - MATHYPRESTRUCT = "hyprestruct" - A matrix type to be used for parallel sparse matrices 294f4bdf6c4SBarry Smith based on the hypre HYPRE_StructMatrix. 295f4bdf6c4SBarry Smith 296f4bdf6c4SBarry Smith Level: intermediate 297f4bdf6c4SBarry Smith 298f4bdf6c4SBarry Smith Notes: Unlike the more general support for blocks in hypre this allows only one block per process and requires the block 299f4bdf6c4SBarry Smith be defined by a DMDA. 300f4bdf6c4SBarry Smith 301950540a4SJed Brown The matrix needs a DMDA associated with it by either a call to MatSetDM() or if the matrix is obtained from DMCreateMatrix() 302f4bdf6c4SBarry Smith 303950540a4SJed Brown .seealso: MatCreate(), PCPFMG, MatSetDM(), DMCreateMatrix() 304f4bdf6c4SBarry Smith M*/ 305f4bdf6c4SBarry Smith 306f4bdf6c4SBarry Smith 307f4bdf6c4SBarry Smith #undef __FUNCT__ 308f4bdf6c4SBarry Smith #define __FUNCT__ "MatSetValuesLocal_HYPREStruct_3d" 309f4bdf6c4SBarry Smith PetscErrorCode MatSetValuesLocal_HYPREStruct_3d(Mat mat,PetscInt nrow,const PetscInt irow[],PetscInt ncol,const PetscInt icol[],const PetscScalar y[],InsertMode addv) 310f4bdf6c4SBarry Smith { 311f4bdf6c4SBarry Smith PetscErrorCode ierr; 312f4bdf6c4SBarry Smith PetscInt i,j,stencil,index[3],row,entries[7]; 313f4bdf6c4SBarry Smith const PetscScalar *values = y; 314f4bdf6c4SBarry Smith Mat_HYPREStruct *ex = (Mat_HYPREStruct*) mat->data; 315f4bdf6c4SBarry Smith 316f4bdf6c4SBarry Smith PetscFunctionBegin; 317f4bdf6c4SBarry Smith for (i=0; i<nrow; i++) { 318f4bdf6c4SBarry Smith for (j=0; j<ncol; j++) { 319f4bdf6c4SBarry Smith stencil = icol[j] - irow[i]; 320f4bdf6c4SBarry Smith if (!stencil) { 321f4bdf6c4SBarry Smith entries[j] = 3; 322f4bdf6c4SBarry Smith } else if (stencil == -1) { 323f4bdf6c4SBarry Smith entries[j] = 2; 324f4bdf6c4SBarry Smith } else if (stencil == 1) { 325f4bdf6c4SBarry Smith entries[j] = 4; 326f4bdf6c4SBarry Smith } else if (stencil == -ex->gnx) { 327f4bdf6c4SBarry Smith entries[j] = 1; 328f4bdf6c4SBarry Smith } else if (stencil == ex->gnx) { 329f4bdf6c4SBarry Smith entries[j] = 5; 330f4bdf6c4SBarry Smith } else if (stencil == -ex->gnxgny) { 331f4bdf6c4SBarry Smith entries[j] = 0; 332f4bdf6c4SBarry Smith } else if (stencil == ex->gnxgny) { 333f4bdf6c4SBarry Smith entries[j] = 6; 334f4bdf6c4SBarry 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); 335f4bdf6c4SBarry Smith } 336f4bdf6c4SBarry Smith row = ex->gindices[irow[i]] - ex->rstart; 337f4bdf6c4SBarry Smith index[0] = ex->xs + (row % ex->nx); 338f4bdf6c4SBarry Smith index[1] = ex->ys + ((row/ex->nx) % ex->ny); 339f4bdf6c4SBarry Smith index[2] = ex->zs + (row/(ex->nxny)); 340f4bdf6c4SBarry Smith if (addv == ADD_VALUES) { 341f4bdf6c4SBarry Smith PetscStackCallHypre(0,HYPRE_StructMatrixAddToValues,(ex->hmat,index,ncol,entries,(PetscScalar*)values)); 342f4bdf6c4SBarry Smith } else { 343f4bdf6c4SBarry Smith PetscStackCallHypre(0,HYPRE_StructMatrixSetValues,(ex->hmat,index,ncol,entries,(PetscScalar*)values)); 344f4bdf6c4SBarry Smith } 345f4bdf6c4SBarry Smith values += ncol; 346f4bdf6c4SBarry Smith } 347f4bdf6c4SBarry Smith PetscFunctionReturn(0); 348f4bdf6c4SBarry Smith } 349f4bdf6c4SBarry Smith 350f4bdf6c4SBarry Smith #undef __FUNCT__ 351f4bdf6c4SBarry Smith #define __FUNCT__ "MatZeroRowsLocal_HYPREStruct_3d" 352f4bdf6c4SBarry Smith PetscErrorCode MatZeroRowsLocal_HYPREStruct_3d(Mat mat,PetscInt nrow,const PetscInt irow[],PetscScalar d,Vec x,Vec b) 353f4bdf6c4SBarry Smith { 354f4bdf6c4SBarry Smith PetscErrorCode ierr; 355f4bdf6c4SBarry Smith PetscInt i,index[3],row,entries[7] = {0,1,2,3,4,5,6}; 356f4bdf6c4SBarry Smith PetscScalar values[7]; 357f4bdf6c4SBarry Smith Mat_HYPREStruct *ex = (Mat_HYPREStruct*) mat->data; 358f4bdf6c4SBarry Smith 359f4bdf6c4SBarry Smith PetscFunctionBegin; 360f4bdf6c4SBarry Smith if (x && b) SETERRQ(((PetscObject)mat)->comm,PETSC_ERR_SUP,"No support"); 361f4bdf6c4SBarry Smith ierr = PetscMemzero(values,7*sizeof(PetscScalar));CHKERRQ(ierr); 362f4bdf6c4SBarry Smith values[3] = d; 363f4bdf6c4SBarry Smith for (i=0; i<nrow; i++) { 364f4bdf6c4SBarry Smith row = ex->gindices[irow[i]] - ex->rstart; 365f4bdf6c4SBarry Smith index[0] = ex->xs + (row % ex->nx); 366f4bdf6c4SBarry Smith index[1] = ex->ys + ((row/ex->nx) % ex->ny); 367f4bdf6c4SBarry Smith index[2] = ex->zs + (row/(ex->nxny)); 368f4bdf6c4SBarry Smith PetscStackCallHypre(0,HYPRE_StructMatrixSetValues,(ex->hmat,index,7,entries,values)); 369f4bdf6c4SBarry Smith } 370f4bdf6c4SBarry Smith PetscStackCallHypre(0,HYPRE_StructMatrixAssemble,(ex->hmat)); 371f4bdf6c4SBarry Smith PetscFunctionReturn(0); 372f4bdf6c4SBarry Smith } 373f4bdf6c4SBarry Smith 374f4bdf6c4SBarry Smith #undef __FUNCT__ 375f4bdf6c4SBarry Smith #define __FUNCT__ "MatZeroEntries_HYPREStruct_3d" 376f4bdf6c4SBarry Smith PetscErrorCode MatZeroEntries_HYPREStruct_3d(Mat mat) 377f4bdf6c4SBarry Smith { 378f4bdf6c4SBarry Smith PetscErrorCode ierr; 379f4bdf6c4SBarry Smith PetscInt indices[7] = {0,1,2,3,4,5,6}; 380f4bdf6c4SBarry Smith Mat_HYPREStruct *ex = (Mat_HYPREStruct*) mat->data; 381f4bdf6c4SBarry Smith 382f4bdf6c4SBarry Smith PetscFunctionBegin; 383f4bdf6c4SBarry Smith /* hypre has no public interface to do this */ 384f4bdf6c4SBarry Smith PetscStackCallHypre(0,hypre_StructMatrixClearBoxValues,(ex->hmat,&ex->hbox,7,indices,0,1)); 385f4bdf6c4SBarry Smith PetscStackCallHypre(0,HYPRE_StructMatrixAssemble,(ex->hmat)); 386f4bdf6c4SBarry Smith PetscFunctionReturn(0); 387f4bdf6c4SBarry Smith } 388f4bdf6c4SBarry Smith 389f4bdf6c4SBarry Smith #undef __FUNCT__ 39095ee5b0eSBarry Smith #define __FUNCT__ "MatSetDM_HYPREStruct" 39195ee5b0eSBarry Smith PetscErrorCode MatSetDM_HYPREStruct(Mat mat,DM da) 392f4bdf6c4SBarry Smith { 393f4bdf6c4SBarry Smith PetscErrorCode ierr; 394f4bdf6c4SBarry Smith Mat_HYPREStruct *ex = (Mat_HYPREStruct*) mat->data; 395f4bdf6c4SBarry Smith PetscInt dim,dof,sw[3],nx,ny,nz; 396f4bdf6c4SBarry Smith int ilower[3],iupper[3],ssize,i; 3977ab44359SJed Brown DMDABoundaryType px,py,pz; 398f4bdf6c4SBarry Smith DMDAStencilType st; 399f4bdf6c4SBarry Smith 400f4bdf6c4SBarry Smith PetscFunctionBegin; 401f4bdf6c4SBarry Smith ex->da = da; 402f4bdf6c4SBarry Smith ierr = PetscObjectReference((PetscObject)da);CHKERRQ(ierr); 403f4bdf6c4SBarry Smith 4047ab44359SJed Brown ierr = DMDAGetInfo(ex->da,&dim,0,0,0,0,0,0,&dof,&sw[0],&px,&py,&pz,&st);CHKERRQ(ierr); 405f4bdf6c4SBarry Smith ierr = DMDAGetCorners(ex->da,&ilower[0],&ilower[1],&ilower[2],&iupper[0],&iupper[1],&iupper[2]);CHKERRQ(ierr); 406f4bdf6c4SBarry Smith iupper[0] += ilower[0] - 1; 407f4bdf6c4SBarry Smith iupper[1] += ilower[1] - 1; 408f4bdf6c4SBarry Smith iupper[2] += ilower[2] - 1; 409f4bdf6c4SBarry Smith 410f4bdf6c4SBarry Smith /* the hypre_Box is used to zero out the matrix entries in MatZeroValues() */ 411f4bdf6c4SBarry Smith ex->hbox.imin[0] = ilower[0]; 412f4bdf6c4SBarry Smith ex->hbox.imin[1] = ilower[1]; 413f4bdf6c4SBarry Smith ex->hbox.imin[2] = ilower[2]; 414f4bdf6c4SBarry Smith ex->hbox.imax[0] = iupper[0]; 415f4bdf6c4SBarry Smith ex->hbox.imax[1] = iupper[1]; 416f4bdf6c4SBarry Smith ex->hbox.imax[2] = iupper[2]; 417f4bdf6c4SBarry Smith 418f4bdf6c4SBarry Smith /* create the hypre grid object and set its information */ 419f4bdf6c4SBarry Smith if (dof > 1) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"Currently only support for scalar problems"); 4207ab44359SJed Brown if (px || py || pz) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"Ask us to add periodic support by calling HYPRE_StructGridSetPeriodic()"); 421f4bdf6c4SBarry Smith PetscStackCallHypre(0,HYPRE_StructGridCreate,(ex->hcomm,dim,&ex->hgrid)); 422f4bdf6c4SBarry Smith PetscStackCallHypre(0,HYPRE_StructGridSetExtents,(ex->hgrid,ilower,iupper)); 423f4bdf6c4SBarry Smith PetscStackCallHypre(0,HYPRE_StructGridAssemble,(ex->hgrid)); 424f4bdf6c4SBarry Smith 425f4bdf6c4SBarry Smith sw[1] = sw[0]; 426f4bdf6c4SBarry Smith sw[2] = sw[1]; 427f4bdf6c4SBarry Smith PetscStackCallHypre(0,HYPRE_StructGridSetNumGhost,(ex->hgrid,sw)); 428f4bdf6c4SBarry Smith 429f4bdf6c4SBarry Smith /* create the hypre stencil object and set its information */ 430f4bdf6c4SBarry Smith if (sw[0] > 1) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"Ask us to add support for wider stencils"); 431f4bdf6c4SBarry Smith if (st == DMDA_STENCIL_BOX) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"Ask us to add support for box stencils"); 432f4bdf6c4SBarry Smith if (dim == 1) { 433f4bdf6c4SBarry Smith int offsets[3][1] = {{-1},{0},{1}}; 434f4bdf6c4SBarry Smith ssize = 3; 435f4bdf6c4SBarry Smith PetscStackCallHypre(0,HYPRE_StructStencilCreate,(dim,ssize,&ex->hstencil)); 436f4bdf6c4SBarry Smith for (i=0; i<ssize; i++) { 437f4bdf6c4SBarry Smith PetscStackCallHypre(0,HYPRE_StructStencilSetElement,(ex->hstencil,i,offsets[i])); 438f4bdf6c4SBarry Smith } 439f4bdf6c4SBarry Smith } else if (dim == 2) { 440f4bdf6c4SBarry Smith int offsets[5][2] = {{0,-1},{-1,0},{0,0},{1,0},{0,1}}; 441f4bdf6c4SBarry Smith ssize = 5; 442f4bdf6c4SBarry Smith PetscStackCallHypre(0,HYPRE_StructStencilCreate,(dim,ssize,&ex->hstencil)); 443f4bdf6c4SBarry Smith for (i=0; i<ssize; i++) { 444f4bdf6c4SBarry Smith PetscStackCallHypre(0,HYPRE_StructStencilSetElement,(ex->hstencil,i,offsets[i])); 445f4bdf6c4SBarry Smith } 446f4bdf6c4SBarry Smith } else if (dim == 3) { 447f4bdf6c4SBarry Smith int offsets[7][3] = {{0,0,-1},{0,-1,0},{-1,0,0},{0,0,0},{1,0,0},{0,1,0},{0,0,1}}; 448f4bdf6c4SBarry Smith ssize = 7; 449f4bdf6c4SBarry Smith PetscStackCallHypre(0,HYPRE_StructStencilCreate,(dim,ssize,&ex->hstencil)); 450f4bdf6c4SBarry Smith for (i=0; i<ssize; i++) { 451f4bdf6c4SBarry Smith PetscStackCallHypre(0,HYPRE_StructStencilSetElement,(ex->hstencil,i,offsets[i])); 452f4bdf6c4SBarry Smith } 453f4bdf6c4SBarry Smith } 454f4bdf6c4SBarry Smith 455f4bdf6c4SBarry Smith /* create the HYPRE vector for rhs and solution */ 456f4bdf6c4SBarry Smith PetscStackCallHypre(0,HYPRE_StructVectorCreate,(ex->hcomm,ex->hgrid,&ex->hb)); 457f4bdf6c4SBarry Smith PetscStackCallHypre(0,HYPRE_StructVectorCreate,(ex->hcomm,ex->hgrid,&ex->hx)); 458f4bdf6c4SBarry Smith PetscStackCallHypre(0,HYPRE_StructVectorInitialize,(ex->hb)); 459f4bdf6c4SBarry Smith PetscStackCallHypre(0,HYPRE_StructVectorInitialize,(ex->hx)); 460f4bdf6c4SBarry Smith PetscStackCallHypre(0,HYPRE_StructVectorAssemble,(ex->hb)); 461f4bdf6c4SBarry Smith PetscStackCallHypre(0,HYPRE_StructVectorAssemble,(ex->hx)); 462f4bdf6c4SBarry Smith 463f4bdf6c4SBarry Smith /* create the hypre matrix object and set its information */ 464f4bdf6c4SBarry Smith PetscStackCallHypre(0,HYPRE_StructMatrixCreate,(ex->hcomm,ex->hgrid,ex->hstencil,&ex->hmat)); 465f4bdf6c4SBarry Smith PetscStackCallHypre(0,HYPRE_StructGridDestroy,(ex->hgrid)); 466f4bdf6c4SBarry Smith PetscStackCallHypre(0,HYPRE_StructStencilDestroy,(ex->hstencil)); 467f4bdf6c4SBarry Smith if (ex->needsinitialization) { 468f4bdf6c4SBarry Smith PetscStackCallHypre(0,HYPRE_StructMatrixInitialize,(ex->hmat)); 469f4bdf6c4SBarry Smith ex->needsinitialization = PETSC_FALSE; 470f4bdf6c4SBarry Smith } 471f4bdf6c4SBarry Smith 472f4bdf6c4SBarry Smith /* set the global and local sizes of the matrix */ 473f4bdf6c4SBarry Smith ierr = DMDAGetCorners(da,0,0,0,&nx,&ny,&nz);CHKERRQ(ierr); 474f4bdf6c4SBarry Smith ierr = MatSetSizes(mat,dof*nx*ny*nz,dof*nx*ny*nz,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr); 475f4bdf6c4SBarry Smith ierr = PetscLayoutSetBlockSize(mat->rmap,1);CHKERRQ(ierr); 476f4bdf6c4SBarry Smith ierr = PetscLayoutSetBlockSize(mat->cmap,1);CHKERRQ(ierr); 477f4bdf6c4SBarry Smith ierr = PetscLayoutSetUp(mat->rmap);CHKERRQ(ierr); 478f4bdf6c4SBarry Smith ierr = PetscLayoutSetUp(mat->cmap);CHKERRQ(ierr); 479f4bdf6c4SBarry Smith 480f4bdf6c4SBarry Smith if (dim == 3) { 481f4bdf6c4SBarry Smith mat->ops->setvalueslocal = MatSetValuesLocal_HYPREStruct_3d; 482f4bdf6c4SBarry Smith mat->ops->zerorowslocal = MatZeroRowsLocal_HYPREStruct_3d; 483f4bdf6c4SBarry Smith mat->ops->zeroentries = MatZeroEntries_HYPREStruct_3d; 484f4bdf6c4SBarry Smith ierr = MatZeroEntries_HYPREStruct_3d(mat);CHKERRQ(ierr); 485f4bdf6c4SBarry Smith } else SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"Only support for 3d DMDA currently"); 486f4bdf6c4SBarry Smith 487f4bdf6c4SBarry Smith /* get values that will be used repeatedly in MatSetValuesLocal() and MatZeroRowsLocal() repeatedly */ 488f4bdf6c4SBarry Smith ierr = MatGetOwnershipRange(mat,&ex->rstart,PETSC_NULL);CHKERRQ(ierr); 489f4bdf6c4SBarry Smith ierr = DMDAGetGlobalIndices(ex->da,PETSC_NULL,&ex->gindices);CHKERRQ(ierr); 490f4bdf6c4SBarry Smith ierr = DMDAGetGhostCorners(ex->da,0,0,0,&ex->gnx,&ex->gnxgny,0);CHKERRQ(ierr); 491f4bdf6c4SBarry Smith ex->gnxgny *= ex->gnx; 492f4bdf6c4SBarry Smith ierr = DMDAGetCorners(ex->da,&ex->xs,&ex->ys,&ex->zs,&ex->nx,&ex->ny,0);CHKERRQ(ierr); 493f4bdf6c4SBarry Smith ex->nxny = ex->nx*ex->ny; 494f4bdf6c4SBarry Smith PetscFunctionReturn(0); 495f4bdf6c4SBarry Smith } 496f4bdf6c4SBarry Smith 497f4bdf6c4SBarry Smith #undef __FUNCT__ 498f4bdf6c4SBarry Smith #define __FUNCT__ "MatMult_HYPREStruct" 499f4bdf6c4SBarry Smith PetscErrorCode MatMult_HYPREStruct(Mat A,Vec x,Vec y) 500f4bdf6c4SBarry Smith { 501f4bdf6c4SBarry Smith PetscErrorCode ierr; 502f4bdf6c4SBarry Smith PetscScalar *xx,*yy; 503f4bdf6c4SBarry Smith int ilower[3],iupper[3]; 504f4bdf6c4SBarry Smith Mat_HYPREStruct *mx = (Mat_HYPREStruct *)(A->data); 505f4bdf6c4SBarry Smith 506f4bdf6c4SBarry Smith PetscFunctionBegin; 507f4bdf6c4SBarry Smith ierr = DMDAGetCorners(mx->da,&ilower[0],&ilower[1],&ilower[2],&iupper[0],&iupper[1],&iupper[2]);CHKERRQ(ierr); 508f4bdf6c4SBarry Smith iupper[0] += ilower[0] - 1; 509f4bdf6c4SBarry Smith iupper[1] += ilower[1] - 1; 510f4bdf6c4SBarry Smith iupper[2] += ilower[2] - 1; 511f4bdf6c4SBarry Smith 512f4bdf6c4SBarry Smith /* copy x values over to hypre */ 513f4bdf6c4SBarry Smith PetscStackCallHypre(0,HYPRE_StructVectorSetConstantValues,(mx->hb,0.0)); 514f4bdf6c4SBarry Smith ierr = VecGetArray(x,&xx);CHKERRQ(ierr); 515f4bdf6c4SBarry Smith PetscStackCallHypre(0,HYPRE_StructVectorSetBoxValues,(mx->hb,ilower,iupper,xx)); 516f4bdf6c4SBarry Smith ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr); 517f4bdf6c4SBarry Smith PetscStackCallHypre(0,HYPRE_StructVectorAssemble,(mx->hb)); 518f4bdf6c4SBarry Smith PetscStackCallHypre(0,HYPRE_StructMatrixMatvec,(1.0,mx->hmat,mx->hb,0.0,mx->hx)); 519f4bdf6c4SBarry Smith 520f4bdf6c4SBarry Smith /* copy solution values back to PETSc */ 521f4bdf6c4SBarry Smith ierr = VecGetArray(y,&yy);CHKERRQ(ierr); 522f4bdf6c4SBarry Smith PetscStackCallHypre(0,HYPRE_StructVectorGetBoxValues,(mx->hx,ilower,iupper,yy)); 523f4bdf6c4SBarry Smith ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr); 524f4bdf6c4SBarry Smith PetscFunctionReturn(0); 525f4bdf6c4SBarry Smith } 526f4bdf6c4SBarry Smith 527f4bdf6c4SBarry Smith #undef __FUNCT__ 528f4bdf6c4SBarry Smith #define __FUNCT__ "MatAssemblyEnd_HYPREStruct" 529f4bdf6c4SBarry Smith PetscErrorCode MatAssemblyEnd_HYPREStruct(Mat mat,MatAssemblyType mode) 530f4bdf6c4SBarry Smith { 531f4bdf6c4SBarry Smith Mat_HYPREStruct *ex = (Mat_HYPREStruct*) mat->data; 532f4bdf6c4SBarry Smith PetscErrorCode ierr; 533f4bdf6c4SBarry Smith 534f4bdf6c4SBarry Smith PetscFunctionBegin; 535f4bdf6c4SBarry Smith PetscStackCallHypre(0,HYPRE_StructMatrixAssemble,(ex->hmat)); 536f4bdf6c4SBarry Smith /* PetscStackCallHypre(0,HYPRE_StructMatrixPrint,("dummy",ex->hmat,0)); */ 537f4bdf6c4SBarry Smith PetscFunctionReturn(0); 538f4bdf6c4SBarry Smith } 539f4bdf6c4SBarry Smith 540f4bdf6c4SBarry Smith #undef __FUNCT__ 541f4bdf6c4SBarry Smith #define __FUNCT__ "MatZeroEntries_HYPREStruct" 542f4bdf6c4SBarry Smith PetscErrorCode MatZeroEntries_HYPREStruct(Mat mat) 543f4bdf6c4SBarry Smith { 544f4bdf6c4SBarry Smith PetscFunctionBegin; 545f4bdf6c4SBarry Smith /* before the DMDA is set to the matrix the zero doesn't need to do anything */ 546f4bdf6c4SBarry Smith PetscFunctionReturn(0); 547f4bdf6c4SBarry Smith } 548f4bdf6c4SBarry Smith 549f4bdf6c4SBarry Smith 550f4bdf6c4SBarry Smith #undef __FUNCT__ 551f4bdf6c4SBarry Smith #define __FUNCT__ "MatDestroy_HYPREStruct" 552f4bdf6c4SBarry Smith PetscErrorCode MatDestroy_HYPREStruct(Mat mat) 553f4bdf6c4SBarry Smith { 554f4bdf6c4SBarry Smith Mat_HYPREStruct *ex = (Mat_HYPREStruct*) mat->data; 555f4bdf6c4SBarry Smith PetscErrorCode ierr; 556f4bdf6c4SBarry Smith 557f4bdf6c4SBarry Smith PetscFunctionBegin; 558f4bdf6c4SBarry Smith PetscStackCallHypre(0,HYPRE_StructMatrixDestroy,(ex->hmat)); 559f4bdf6c4SBarry Smith PetscStackCallHypre(0,HYPRE_StructVectorDestroy,(ex->hx)); 560f4bdf6c4SBarry Smith PetscStackCallHypre(0,HYPRE_StructVectorDestroy,(ex->hb)); 561f4bdf6c4SBarry Smith PetscFunctionReturn(0); 562f4bdf6c4SBarry Smith } 563f4bdf6c4SBarry Smith 564f4bdf6c4SBarry Smith 565f4bdf6c4SBarry Smith EXTERN_C_BEGIN 566f4bdf6c4SBarry Smith #undef __FUNCT__ 567f4bdf6c4SBarry Smith #define __FUNCT__ "MatCreate_HYPREStruct" 568f4bdf6c4SBarry Smith PetscErrorCode MatCreate_HYPREStruct(Mat B) 569f4bdf6c4SBarry Smith { 570f4bdf6c4SBarry Smith Mat_HYPREStruct *ex; 571f4bdf6c4SBarry Smith PetscErrorCode ierr; 572f4bdf6c4SBarry Smith 573f4bdf6c4SBarry Smith PetscFunctionBegin; 574f4bdf6c4SBarry Smith ierr = PetscNewLog(B,Mat_HYPREStruct,&ex);CHKERRQ(ierr); 575f4bdf6c4SBarry Smith B->data = (void*)ex; 576f4bdf6c4SBarry Smith B->rmap->bs = 1; 577f4bdf6c4SBarry Smith B->assembled = PETSC_FALSE; 578f4bdf6c4SBarry Smith 579f4bdf6c4SBarry Smith B->insertmode = NOT_SET_VALUES; 580f4bdf6c4SBarry Smith 581f4bdf6c4SBarry Smith B->ops->assemblyend = MatAssemblyEnd_HYPREStruct; 582f4bdf6c4SBarry Smith B->ops->mult = MatMult_HYPREStruct; 583f4bdf6c4SBarry Smith B->ops->zeroentries = MatZeroEntries_HYPREStruct; 584f4bdf6c4SBarry Smith B->ops->destroy = MatDestroy_HYPREStruct; 585f4bdf6c4SBarry Smith 586f4bdf6c4SBarry Smith ex->needsinitialization = PETSC_TRUE; 587f4bdf6c4SBarry Smith 588f4bdf6c4SBarry Smith ierr = MPI_Comm_dup(((PetscObject)B)->comm,&(ex->hcomm));CHKERRQ(ierr); 58995ee5b0eSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSetDM_C","MatSetDM_HYPREStruct",MatSetDM_HYPREStruct);CHKERRQ(ierr); 590f4bdf6c4SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)B,MATHYPRESTRUCT);CHKERRQ(ierr); 591f4bdf6c4SBarry Smith PetscFunctionReturn(0); 592f4bdf6c4SBarry Smith } 593f4bdf6c4SBarry Smith EXTERN_C_END 594f4bdf6c4SBarry Smith 595f4bdf6c4SBarry Smith 596f4bdf6c4SBarry Smith /*MC 597f4bdf6c4SBarry Smith MATHYPRESSTRUCT - MATHYPRESSTRUCT = "hypresstruct" - A matrix type to be used for parallel sparse matrices 598f4bdf6c4SBarry Smith based on the hypre HYPRE_SStructMatrix. 599f4bdf6c4SBarry Smith 600f4bdf6c4SBarry Smith 601f4bdf6c4SBarry Smith Level: intermediate 602f4bdf6c4SBarry Smith 603f4bdf6c4SBarry Smith Notes: Unlike hypre's general semi-struct object consisting of a collection of structured-grid objects and unstructured 604f4bdf6c4SBarry Smith grid objects, we will restrict the semi-struct objects to consist of only structured-grid components. 605f4bdf6c4SBarry Smith 606f4bdf6c4SBarry 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 607f4bdf6c4SBarry Smith be defined by a DMDA. 608f4bdf6c4SBarry Smith 609950540a4SJed Brown The matrix needs a DMDA associated with it by either a call to MatSetDM() or if the matrix is obtained from DMCreateMatrix() 610f4bdf6c4SBarry Smith 611f4bdf6c4SBarry Smith M*/ 612f4bdf6c4SBarry Smith 613f4bdf6c4SBarry Smith #undef __FUNCT__ 614f4bdf6c4SBarry Smith #define __FUNCT__ "MatSetValuesLocal_HYPRESStruct_3d" 615f4bdf6c4SBarry Smith PetscErrorCode MatSetValuesLocal_HYPRESStruct_3d(Mat mat,PetscInt nrow,const PetscInt irow[],PetscInt ncol,const PetscInt icol[],const PetscScalar y[],InsertMode addv) 616f4bdf6c4SBarry Smith { 617f4bdf6c4SBarry Smith PetscErrorCode ierr; 618f4bdf6c4SBarry Smith PetscInt i,j,stencil,index[3]; 619f4bdf6c4SBarry Smith const PetscScalar *values = y; 620f4bdf6c4SBarry Smith Mat_HYPRESStruct *ex = (Mat_HYPRESStruct*) mat->data; 621f4bdf6c4SBarry Smith 622f4bdf6c4SBarry Smith int part= 0; /* Petsc sstruct interface only allows 1 part */ 623f4bdf6c4SBarry Smith int ordering; 624f4bdf6c4SBarry Smith int grid_rank, to_grid_rank; 625f4bdf6c4SBarry Smith int var_type, to_var_type; 626f4bdf6c4SBarry Smith int to_var_entry = 0; 627f4bdf6c4SBarry Smith 628f4bdf6c4SBarry Smith int nvars= ex->nvars; 629f4bdf6c4SBarry Smith PetscInt row,*entries; 630f4bdf6c4SBarry Smith 631f4bdf6c4SBarry Smith PetscFunctionBegin; 632f4bdf6c4SBarry Smith ierr = PetscMalloc(7*nvars*sizeof(PetscInt),&entries);CHKERRQ(ierr); 633f4bdf6c4SBarry Smith ordering= ex-> dofs_order; /* ordering= 0 nodal ordering 634f4bdf6c4SBarry Smith 1 variable ordering */ 635f4bdf6c4SBarry Smith /* stencil entries are orderer by variables: var0_stencil0, var0_stencil1, ..., var0_stencil6, var1_stencil0, var1_stencil1, ... */ 636f4bdf6c4SBarry Smith 637f4bdf6c4SBarry Smith if (!ordering) { 638f4bdf6c4SBarry Smith for (i=0; i<nrow; i++) { 639f4bdf6c4SBarry Smith grid_rank= irow[i]/nvars; 640f4bdf6c4SBarry Smith var_type = (irow[i] % nvars); 641f4bdf6c4SBarry Smith 642f4bdf6c4SBarry Smith for (j=0; j<ncol; j++) { 643f4bdf6c4SBarry Smith to_grid_rank= icol[j]/nvars; 644f4bdf6c4SBarry Smith to_var_type = (icol[j] % nvars); 645f4bdf6c4SBarry Smith 646f4bdf6c4SBarry Smith to_var_entry= to_var_entry*7; 647f4bdf6c4SBarry Smith entries[j]= to_var_entry; 648f4bdf6c4SBarry Smith 649f4bdf6c4SBarry Smith stencil = to_grid_rank-grid_rank; 650f4bdf6c4SBarry Smith if (!stencil) { 651f4bdf6c4SBarry Smith entries[j] += 3; 652f4bdf6c4SBarry Smith } else if (stencil == -1) { 653f4bdf6c4SBarry Smith entries[j] += 2; 654f4bdf6c4SBarry Smith } else if (stencil == 1) { 655f4bdf6c4SBarry Smith entries[j] += 4; 656f4bdf6c4SBarry Smith } else if (stencil == -ex->gnx) { 657f4bdf6c4SBarry Smith entries[j] += 1; 658f4bdf6c4SBarry Smith } else if (stencil == ex->gnx) { 659f4bdf6c4SBarry Smith entries[j] += 5; 660f4bdf6c4SBarry Smith } else if (stencil == -ex->gnxgny) { 661f4bdf6c4SBarry Smith entries[j] += 0; 662f4bdf6c4SBarry Smith } else if (stencil == ex->gnxgny) { 663f4bdf6c4SBarry Smith entries[j] += 6; 664f4bdf6c4SBarry 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); 665f4bdf6c4SBarry Smith } 666f4bdf6c4SBarry Smith 667f4bdf6c4SBarry Smith row = ex->gindices[grid_rank] - ex->rstart; 668f4bdf6c4SBarry Smith index[0] = ex->xs + (row % ex->nx); 669f4bdf6c4SBarry Smith index[1] = ex->ys + ((row/ex->nx) % ex->ny); 670f4bdf6c4SBarry Smith index[2] = ex->zs + (row/(ex->nxny)); 671f4bdf6c4SBarry Smith 672f4bdf6c4SBarry Smith if (addv == ADD_VALUES) { 673f4bdf6c4SBarry Smith PetscStackCallHypre(0,HYPRE_SStructMatrixAddToValues,(ex->ss_mat,part,index,var_type,ncol,entries,(PetscScalar*)values)); 674f4bdf6c4SBarry Smith } else { 675f4bdf6c4SBarry Smith PetscStackCallHypre(0,HYPRE_SStructMatrixSetValues,(ex->ss_mat,part,index,var_type,ncol,entries,(PetscScalar*)values)); 676f4bdf6c4SBarry Smith } 677f4bdf6c4SBarry Smith values += ncol; 678f4bdf6c4SBarry Smith } 679f4bdf6c4SBarry Smith } else { 680f4bdf6c4SBarry Smith for (i=0; i<nrow; i++) { 681f4bdf6c4SBarry Smith var_type = irow[i]/(ex->gnxgnygnz); 682f4bdf6c4SBarry Smith grid_rank= irow[i] - var_type*(ex->gnxgnygnz); 683f4bdf6c4SBarry Smith 684f4bdf6c4SBarry Smith for (j=0; j<ncol; j++) { 685f4bdf6c4SBarry Smith to_var_type = icol[j]/(ex->gnxgnygnz); 686f4bdf6c4SBarry Smith to_grid_rank= icol[j] - to_var_type*(ex->gnxgnygnz); 687f4bdf6c4SBarry Smith 688f4bdf6c4SBarry Smith to_var_entry= to_var_entry*7; 689f4bdf6c4SBarry Smith entries[j]= to_var_entry; 690f4bdf6c4SBarry Smith 691f4bdf6c4SBarry Smith stencil = to_grid_rank-grid_rank; 692f4bdf6c4SBarry Smith if (!stencil) { 693f4bdf6c4SBarry Smith entries[j] += 3; 694f4bdf6c4SBarry Smith } else if (stencil == -1) { 695f4bdf6c4SBarry Smith entries[j] += 2; 696f4bdf6c4SBarry Smith } else if (stencil == 1) { 697f4bdf6c4SBarry Smith entries[j] += 4; 698f4bdf6c4SBarry Smith } else if (stencil == -ex->gnx) { 699f4bdf6c4SBarry Smith entries[j] += 1; 700f4bdf6c4SBarry Smith } else if (stencil == ex->gnx) { 701f4bdf6c4SBarry Smith entries[j] += 5; 702f4bdf6c4SBarry Smith } else if (stencil == -ex->gnxgny) { 703f4bdf6c4SBarry Smith entries[j] += 0; 704f4bdf6c4SBarry Smith } else if (stencil == ex->gnxgny) { 705f4bdf6c4SBarry Smith entries[j] += 6; 706f4bdf6c4SBarry 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); 707f4bdf6c4SBarry Smith } 708f4bdf6c4SBarry Smith 709f4bdf6c4SBarry Smith row = ex->gindices[grid_rank] - ex->rstart; 710f4bdf6c4SBarry Smith index[0] = ex->xs + (row % ex->nx); 711f4bdf6c4SBarry Smith index[1] = ex->ys + ((row/ex->nx) % ex->ny); 712f4bdf6c4SBarry Smith index[2] = ex->zs + (row/(ex->nxny)); 713f4bdf6c4SBarry Smith 714f4bdf6c4SBarry Smith if (addv == ADD_VALUES) { 715f4bdf6c4SBarry Smith PetscStackCallHypre(0,HYPRE_SStructMatrixAddToValues,(ex->ss_mat,part,index,var_type,ncol,entries,(PetscScalar*)values)); 716f4bdf6c4SBarry Smith } else { 717f4bdf6c4SBarry Smith PetscStackCallHypre(0,HYPRE_SStructMatrixSetValues,(ex->ss_mat,part,index,var_type,ncol,entries,(PetscScalar*)values)); 718f4bdf6c4SBarry Smith } 719f4bdf6c4SBarry Smith values += ncol; 720f4bdf6c4SBarry Smith } 721f4bdf6c4SBarry Smith } 722f4bdf6c4SBarry Smith ierr = PetscFree(entries);CHKERRQ(ierr); 723f4bdf6c4SBarry Smith PetscFunctionReturn(0); 724f4bdf6c4SBarry Smith } 725f4bdf6c4SBarry Smith 726f4bdf6c4SBarry Smith #undef __FUNCT__ 727f4bdf6c4SBarry Smith #define __FUNCT__ "MatZeroRowsLocal_HYPRESStruct_3d" 728f4bdf6c4SBarry Smith PetscErrorCode MatZeroRowsLocal_HYPRESStruct_3d(Mat mat,PetscInt nrow,const PetscInt irow[],PetscScalar d,Vec x,Vec b) 729f4bdf6c4SBarry Smith { 730f4bdf6c4SBarry Smith PetscErrorCode ierr; 731f4bdf6c4SBarry Smith PetscInt i,index[3]; 732f4bdf6c4SBarry Smith PetscScalar **values; 733f4bdf6c4SBarry Smith Mat_HYPRESStruct *ex = (Mat_HYPRESStruct*) mat->data; 734f4bdf6c4SBarry Smith 735f4bdf6c4SBarry Smith int part= 0; /* Petsc sstruct interface only allows 1 part */ 736f4bdf6c4SBarry Smith int ordering= ex->dofs_order; 737f4bdf6c4SBarry Smith int grid_rank; 738f4bdf6c4SBarry Smith int var_type; 739f4bdf6c4SBarry Smith int nvars= ex->nvars; 740f4bdf6c4SBarry Smith PetscInt row,*entries; 741f4bdf6c4SBarry Smith 742f4bdf6c4SBarry Smith PetscFunctionBegin; 743f4bdf6c4SBarry Smith if (x && b) SETERRQ(((PetscObject)mat)->comm,PETSC_ERR_SUP,"No support"); 744f4bdf6c4SBarry Smith ierr = PetscMalloc(7*nvars*sizeof(PetscInt),&entries);CHKERRQ(ierr); 745f4bdf6c4SBarry Smith 746f4bdf6c4SBarry Smith ierr = PetscMalloc(nvars*sizeof(PetscScalar *),&values);CHKERRQ(ierr); 747f4bdf6c4SBarry Smith ierr = PetscMalloc(7*nvars*nvars*sizeof(PetscScalar),&values[0]);CHKERRQ(ierr); 748f4bdf6c4SBarry Smith for (i=1; i<nvars; i++) { 749f4bdf6c4SBarry Smith values[i] = values[i-1] + nvars*7; 750f4bdf6c4SBarry Smith } 751f4bdf6c4SBarry Smith 752f4bdf6c4SBarry Smith for (i=0; i< nvars; i++) { 753f4bdf6c4SBarry Smith ierr = PetscMemzero(values[i],nvars*7*sizeof(PetscScalar));CHKERRQ(ierr); 754f4bdf6c4SBarry Smith *(values[i]+3)= d; 755f4bdf6c4SBarry Smith } 756f4bdf6c4SBarry Smith 757f4bdf6c4SBarry Smith for (i= 0; i< nvars*7; i++) { 758f4bdf6c4SBarry Smith entries[i]= i; 759f4bdf6c4SBarry Smith } 760f4bdf6c4SBarry Smith 761f4bdf6c4SBarry Smith if (!ordering) { 762f4bdf6c4SBarry Smith for (i=0; i<nrow; i++) { 763f4bdf6c4SBarry Smith grid_rank= irow[i]/nvars; 764f4bdf6c4SBarry Smith var_type = (irow[i] % nvars); 765f4bdf6c4SBarry Smith 766f4bdf6c4SBarry Smith row = ex->gindices[grid_rank] - ex->rstart; 767f4bdf6c4SBarry Smith index[0] = ex->xs + (row % ex->nx); 768f4bdf6c4SBarry Smith index[1] = ex->ys + ((row/ex->nx) % ex->ny); 769f4bdf6c4SBarry Smith index[2] = ex->zs + (row/(ex->nxny)); 770f4bdf6c4SBarry Smith PetscStackCallHypre(0,HYPRE_SStructMatrixSetValues,(ex->ss_mat,part,index,var_type,7*nvars,entries,values[var_type])); 771f4bdf6c4SBarry Smith } 772f4bdf6c4SBarry Smith } else { 773f4bdf6c4SBarry Smith for (i=0; i<nrow; i++) { 774f4bdf6c4SBarry Smith var_type = irow[i]/(ex->gnxgnygnz); 775f4bdf6c4SBarry Smith grid_rank= irow[i] - var_type*(ex->gnxgnygnz); 776f4bdf6c4SBarry Smith 777f4bdf6c4SBarry Smith row = ex->gindices[grid_rank] - ex->rstart; 778f4bdf6c4SBarry Smith index[0] = ex->xs + (row % ex->nx); 779f4bdf6c4SBarry Smith index[1] = ex->ys + ((row/ex->nx) % ex->ny); 780f4bdf6c4SBarry Smith index[2] = ex->zs + (row/(ex->nxny)); 781f4bdf6c4SBarry Smith PetscStackCallHypre(0,HYPRE_SStructMatrixSetValues,(ex->ss_mat,part,index,var_type,7*nvars,entries,values[var_type])); 782f4bdf6c4SBarry Smith } 783f4bdf6c4SBarry Smith } 784f4bdf6c4SBarry Smith PetscStackCallHypre(0,HYPRE_SStructMatrixAssemble,(ex->ss_mat)); 785f4bdf6c4SBarry Smith ierr = PetscFree(values[0]);CHKERRQ(ierr); 786f4bdf6c4SBarry Smith ierr = PetscFree(values);CHKERRQ(ierr); 787f4bdf6c4SBarry Smith ierr = PetscFree(entries);CHKERRQ(ierr); 788f4bdf6c4SBarry Smith PetscFunctionReturn(0); 789f4bdf6c4SBarry Smith } 790f4bdf6c4SBarry Smith 791f4bdf6c4SBarry Smith #undef __FUNCT__ 792f4bdf6c4SBarry Smith #define __FUNCT__ "MatZeroEntries_HYPRESStruct_3d" 793f4bdf6c4SBarry Smith PetscErrorCode MatZeroEntries_HYPRESStruct_3d(Mat mat) 794f4bdf6c4SBarry Smith { 795f4bdf6c4SBarry Smith PetscErrorCode ierr; 796f4bdf6c4SBarry Smith Mat_HYPRESStruct *ex = (Mat_HYPRESStruct*) mat->data; 797f4bdf6c4SBarry Smith int nvars= ex->nvars; 798f4bdf6c4SBarry Smith int size; 799f4bdf6c4SBarry Smith int part= 0; /* only one part */ 800f4bdf6c4SBarry Smith 801f4bdf6c4SBarry Smith PetscFunctionBegin; 802f4bdf6c4SBarry 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); 803f4bdf6c4SBarry Smith { 804f4bdf6c4SBarry Smith PetscInt i,*entries; 805f4bdf6c4SBarry Smith PetscScalar *values; 806f4bdf6c4SBarry Smith int iupper[3], ilower[3]; 807f4bdf6c4SBarry Smith 808f4bdf6c4SBarry Smith for (i= 0; i< 3; i++) { 809f4bdf6c4SBarry Smith ilower[i]= ex->hbox.imin[i]; 810f4bdf6c4SBarry Smith iupper[i]= ex->hbox.imax[i]; 811f4bdf6c4SBarry Smith } 812f4bdf6c4SBarry Smith 813f4bdf6c4SBarry Smith ierr = PetscMalloc2(nvars*7,PetscInt,&entries,nvars*7*size,PetscScalar,&values);CHKERRQ(ierr); 814f4bdf6c4SBarry Smith for (i= 0; i< nvars*7; i++) { 815f4bdf6c4SBarry Smith entries[i]= i; 816f4bdf6c4SBarry Smith } 817f4bdf6c4SBarry Smith ierr = PetscMemzero(values,nvars*7*size*sizeof(PetscScalar));CHKERRQ(ierr); 818f4bdf6c4SBarry Smith 819f4bdf6c4SBarry Smith for (i= 0; i< nvars; i++) { 820f4bdf6c4SBarry Smith PetscStackCallHypre(0,HYPRE_SStructMatrixSetBoxValues,(ex->ss_mat,part,ilower,iupper,i,nvars*7,entries,values)); 821f4bdf6c4SBarry Smith } 822f4bdf6c4SBarry Smith ierr = PetscFree2(entries,values);CHKERRQ(ierr); 823f4bdf6c4SBarry Smith } 824f4bdf6c4SBarry Smith PetscStackCallHypre(0,HYPRE_SStructMatrixAssemble,(ex->ss_mat)); 825f4bdf6c4SBarry Smith PetscFunctionReturn(0); 826f4bdf6c4SBarry Smith } 827f4bdf6c4SBarry Smith 828f4bdf6c4SBarry Smith 829f4bdf6c4SBarry Smith #undef __FUNCT__ 83095ee5b0eSBarry Smith #define __FUNCT__ "MatSetDM_HYPRESStruct" 83195ee5b0eSBarry Smith PetscErrorCode MatSetDM_HYPRESStruct(Mat mat,DM da) 832f4bdf6c4SBarry Smith { 833f4bdf6c4SBarry Smith PetscErrorCode ierr; 834f4bdf6c4SBarry Smith Mat_HYPRESStruct *ex = (Mat_HYPRESStruct*) mat->data; 835f4bdf6c4SBarry Smith PetscInt dim,dof,sw[3],nx,ny,nz; 836f4bdf6c4SBarry Smith int ilower[3],iupper[3],ssize,i; 8377ab44359SJed Brown DMDABoundaryType px,py,pz; 838f4bdf6c4SBarry Smith DMDAStencilType st; 839f4bdf6c4SBarry Smith int nparts= 1; /* assuming only one part */ 840f4bdf6c4SBarry Smith int part = 0; 841f4bdf6c4SBarry Smith 842f4bdf6c4SBarry Smith PetscFunctionBegin; 843f4bdf6c4SBarry Smith ex->da = da; 844f4bdf6c4SBarry Smith ierr = PetscObjectReference((PetscObject)da);CHKERRQ(ierr); 845f4bdf6c4SBarry Smith 8467ab44359SJed Brown ierr = DMDAGetInfo(ex->da,&dim,0,0,0,0,0,0,&dof,&sw[0],&px,&py,&pz,&st);CHKERRQ(ierr); 847f4bdf6c4SBarry Smith ierr = DMDAGetCorners(ex->da,&ilower[0],&ilower[1],&ilower[2],&iupper[0],&iupper[1],&iupper[2]);CHKERRQ(ierr); 848f4bdf6c4SBarry Smith iupper[0] += ilower[0] - 1; 849f4bdf6c4SBarry Smith iupper[1] += ilower[1] - 1; 850f4bdf6c4SBarry Smith iupper[2] += ilower[2] - 1; 851f4bdf6c4SBarry Smith /* the hypre_Box is used to zero out the matrix entries in MatZeroValues() */ 852f4bdf6c4SBarry Smith ex->hbox.imin[0] = ilower[0]; 853f4bdf6c4SBarry Smith ex->hbox.imin[1] = ilower[1]; 854f4bdf6c4SBarry Smith ex->hbox.imin[2] = ilower[2]; 855f4bdf6c4SBarry Smith ex->hbox.imax[0] = iupper[0]; 856f4bdf6c4SBarry Smith ex->hbox.imax[1] = iupper[1]; 857f4bdf6c4SBarry Smith ex->hbox.imax[2] = iupper[2]; 858f4bdf6c4SBarry Smith 859f4bdf6c4SBarry Smith ex->dofs_order = 0; 860f4bdf6c4SBarry Smith 861f4bdf6c4SBarry Smith /* assuming that the same number of dofs on each gridpoint. Also assume all cell-centred based */ 862f4bdf6c4SBarry Smith ex->nvars= dof; 863f4bdf6c4SBarry Smith 864f4bdf6c4SBarry Smith /* create the hypre grid object and set its information */ 8657ab44359SJed Brown if (px || py || pz) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"Ask us to add periodic support by calling HYPRE_SStructGridSetPeriodic()"); 866f4bdf6c4SBarry Smith PetscStackCallHypre(0,HYPRE_SStructGridCreate,(ex->hcomm,dim,nparts,&ex->ss_grid)); 867f4bdf6c4SBarry Smith 868f4bdf6c4SBarry Smith PetscStackCallHypre(0,HYPRE_SStructGridSetExtents,(ex->ss_grid,part,ex->hbox.imin,ex->hbox.imax)); 869f4bdf6c4SBarry Smith 870f4bdf6c4SBarry Smith { 871f4bdf6c4SBarry Smith HYPRE_SStructVariable *vartypes; 872f4bdf6c4SBarry Smith ierr = PetscMalloc(ex->nvars*sizeof(HYPRE_SStructVariable),&vartypes);CHKERRQ(ierr); 873f4bdf6c4SBarry Smith for (i= 0; i< ex->nvars; i++) { 874f4bdf6c4SBarry Smith vartypes[i]= HYPRE_SSTRUCT_VARIABLE_CELL; 875f4bdf6c4SBarry Smith } 876f4bdf6c4SBarry Smith PetscStackCallHypre(0,HYPRE_SStructGridSetVariables,(ex->ss_grid, part, ex->nvars,vartypes)); 877f4bdf6c4SBarry Smith ierr = PetscFree(vartypes);CHKERRQ(ierr); 878f4bdf6c4SBarry Smith } 879f4bdf6c4SBarry Smith PetscStackCallHypre(0,HYPRE_SStructGridAssemble,(ex->ss_grid)); 880f4bdf6c4SBarry Smith 881f4bdf6c4SBarry Smith sw[1] = sw[0]; 882f4bdf6c4SBarry Smith sw[2] = sw[1]; 883f4bdf6c4SBarry Smith /* PetscStackCallHypre(0,HYPRE_SStructGridSetNumGhost,(ex->ss_grid,sw)); */ 884f4bdf6c4SBarry Smith 885f4bdf6c4SBarry Smith /* create the hypre stencil object and set its information */ 886f4bdf6c4SBarry Smith if (sw[0] > 1) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"Ask us to add support for wider stencils"); 887f4bdf6c4SBarry Smith if (st == DMDA_STENCIL_BOX) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"Ask us to add support for box stencils"); 888f4bdf6c4SBarry Smith 889f4bdf6c4SBarry Smith if (dim == 1) { 890f4bdf6c4SBarry Smith int offsets[3][1] = {{-1},{0},{1}}; 891f4bdf6c4SBarry Smith int j, cnt; 892f4bdf6c4SBarry Smith 893f4bdf6c4SBarry Smith ssize = 3*(ex->nvars); 894f4bdf6c4SBarry Smith PetscStackCallHypre(0,HYPRE_SStructStencilCreate,(dim,ssize,&ex->ss_stencil)); 895f4bdf6c4SBarry Smith cnt= 0; 896f4bdf6c4SBarry Smith for (i= 0; i< (ex->nvars); i++) { 897f4bdf6c4SBarry Smith for (j= 0; j< 3; j++) { 898f4bdf6c4SBarry Smith PetscStackCallHypre(0,HYPRE_SStructStencilSetEntry,(ex->ss_stencil, cnt, offsets[j], i)); 899f4bdf6c4SBarry Smith cnt++; 900f4bdf6c4SBarry Smith } 901f4bdf6c4SBarry Smith } 902f4bdf6c4SBarry Smith } else if (dim == 2) { 903f4bdf6c4SBarry Smith int offsets[5][2] = {{0,-1},{-1,0},{0,0},{1,0},{0,1}}; 904f4bdf6c4SBarry Smith int j, cnt; 905f4bdf6c4SBarry Smith 906f4bdf6c4SBarry Smith ssize = 5*(ex->nvars); 907f4bdf6c4SBarry Smith PetscStackCallHypre(0,HYPRE_SStructStencilCreate,(dim,ssize,&ex->ss_stencil)); 908f4bdf6c4SBarry Smith cnt= 0; 909f4bdf6c4SBarry Smith for (i= 0; i< (ex->nvars); i++) { 910f4bdf6c4SBarry Smith for (j= 0; j< 5; j++) { 911f4bdf6c4SBarry Smith PetscStackCallHypre(0,HYPRE_SStructStencilSetEntry,(ex->ss_stencil, cnt, offsets[j], i)); 912f4bdf6c4SBarry Smith cnt++; 913f4bdf6c4SBarry Smith } 914f4bdf6c4SBarry Smith } 915f4bdf6c4SBarry Smith } else if (dim == 3) { 916f4bdf6c4SBarry Smith int offsets[7][3] = {{0,0,-1},{0,-1,0},{-1,0,0},{0,0,0},{1,0,0},{0,1,0},{0,0,1}}; 917f4bdf6c4SBarry Smith int j, cnt; 918f4bdf6c4SBarry Smith 919f4bdf6c4SBarry Smith ssize = 7*(ex->nvars); 920f4bdf6c4SBarry Smith PetscStackCallHypre(0,HYPRE_SStructStencilCreate,(dim,ssize,&ex->ss_stencil)); 921f4bdf6c4SBarry Smith cnt= 0; 922f4bdf6c4SBarry Smith for (i= 0; i< (ex->nvars); i++) { 923f4bdf6c4SBarry Smith for (j= 0; j< 7; j++) { 924f4bdf6c4SBarry Smith PetscStackCallHypre(0,HYPRE_SStructStencilSetEntry,(ex->ss_stencil, cnt, offsets[j], i)); 925f4bdf6c4SBarry Smith cnt++; 926f4bdf6c4SBarry Smith } 927f4bdf6c4SBarry Smith } 928f4bdf6c4SBarry Smith } 929f4bdf6c4SBarry Smith 930f4bdf6c4SBarry Smith /* create the HYPRE graph */ 931f4bdf6c4SBarry Smith PetscStackCallHypre(0,HYPRE_SStructGraphCreate,(ex->hcomm, ex->ss_grid, &(ex->ss_graph))); 932f4bdf6c4SBarry Smith 933f4bdf6c4SBarry Smith /* set the stencil graph. Note that each variable has the same graph. This means that each 934f4bdf6c4SBarry Smith variable couples to all the other variable and with the same stencil pattern. */ 935f4bdf6c4SBarry Smith for (i= 0; i< (ex->nvars); i++) { 936f4bdf6c4SBarry Smith PetscStackCallHypre(0,HYPRE_SStructGraphSetStencil,(ex->ss_graph,part,i,ex->ss_stencil)); 937f4bdf6c4SBarry Smith } 938f4bdf6c4SBarry Smith PetscStackCallHypre(0,HYPRE_SStructGraphAssemble,(ex->ss_graph)); 939f4bdf6c4SBarry Smith 940f4bdf6c4SBarry Smith /* create the HYPRE sstruct vectors for rhs and solution */ 941f4bdf6c4SBarry Smith PetscStackCallHypre(0,HYPRE_SStructVectorCreate,(ex->hcomm,ex->ss_grid,&ex->ss_b)); 942f4bdf6c4SBarry Smith PetscStackCallHypre(0,HYPRE_SStructVectorCreate,(ex->hcomm,ex->ss_grid,&ex->ss_x)); 943f4bdf6c4SBarry Smith PetscStackCallHypre(0,HYPRE_SStructVectorInitialize,(ex->ss_b)); 944f4bdf6c4SBarry Smith PetscStackCallHypre(0,HYPRE_SStructVectorInitialize,(ex->ss_x)); 945f4bdf6c4SBarry Smith PetscStackCallHypre(0,HYPRE_SStructVectorAssemble,(ex->ss_b)); 946f4bdf6c4SBarry Smith PetscStackCallHypre(0,HYPRE_SStructVectorAssemble,(ex->ss_x)); 947f4bdf6c4SBarry Smith 948f4bdf6c4SBarry Smith /* create the hypre matrix object and set its information */ 949f4bdf6c4SBarry Smith PetscStackCallHypre(0,HYPRE_SStructMatrixCreate,(ex->hcomm,ex->ss_graph,&ex->ss_mat)); 950f4bdf6c4SBarry Smith PetscStackCallHypre(0,HYPRE_SStructGridDestroy,(ex->ss_grid)); 951f4bdf6c4SBarry Smith PetscStackCallHypre(0,HYPRE_SStructStencilDestroy,(ex->ss_stencil)); 952f4bdf6c4SBarry Smith if (ex->needsinitialization) { 953f4bdf6c4SBarry Smith PetscStackCallHypre(0,HYPRE_SStructMatrixInitialize,(ex->ss_mat)); 954f4bdf6c4SBarry Smith ex->needsinitialization = PETSC_FALSE; 955f4bdf6c4SBarry Smith } 956f4bdf6c4SBarry Smith 957f4bdf6c4SBarry Smith /* set the global and local sizes of the matrix */ 958f4bdf6c4SBarry Smith ierr = DMDAGetCorners(da,0,0,0,&nx,&ny,&nz);CHKERRQ(ierr); 959f4bdf6c4SBarry Smith ierr = MatSetSizes(mat,dof*nx*ny*nz,dof*nx*ny*nz,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr); 960f4bdf6c4SBarry Smith ierr = PetscLayoutSetBlockSize(mat->rmap,1);CHKERRQ(ierr); 961f4bdf6c4SBarry Smith ierr = PetscLayoutSetBlockSize(mat->cmap,1);CHKERRQ(ierr); 962f4bdf6c4SBarry Smith ierr = PetscLayoutSetUp(mat->rmap);CHKERRQ(ierr); 963f4bdf6c4SBarry Smith ierr = PetscLayoutSetUp(mat->cmap);CHKERRQ(ierr); 964f4bdf6c4SBarry Smith 965f4bdf6c4SBarry Smith if (dim == 3) { 966f4bdf6c4SBarry Smith mat->ops->setvalueslocal = MatSetValuesLocal_HYPRESStruct_3d; 967f4bdf6c4SBarry Smith mat->ops->zerorowslocal = MatZeroRowsLocal_HYPRESStruct_3d; 968f4bdf6c4SBarry Smith mat->ops->zeroentries = MatZeroEntries_HYPRESStruct_3d; 969f4bdf6c4SBarry Smith ierr = MatZeroEntries_HYPRESStruct_3d(mat);CHKERRQ(ierr); 970f4bdf6c4SBarry Smith } else SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"Only support for 3d DMDA currently"); 971f4bdf6c4SBarry Smith 972f4bdf6c4SBarry Smith /* get values that will be used repeatedly in MatSetValuesLocal() and MatZeroRowsLocal() repeatedly */ 973f4bdf6c4SBarry Smith ierr = MatGetOwnershipRange(mat,&ex->rstart,PETSC_NULL);CHKERRQ(ierr); 974f4bdf6c4SBarry Smith ierr = DMDAGetGlobalIndices(ex->da,PETSC_NULL,&ex->gindices);CHKERRQ(ierr); 975f4bdf6c4SBarry Smith ierr = DMDAGetGhostCorners(ex->da,0,0,0,&ex->gnx,&ex->gnxgny,&ex->gnxgnygnz);CHKERRQ(ierr); 976f4bdf6c4SBarry Smith ex->gnxgny *= ex->gnx; 977f4bdf6c4SBarry Smith ex->gnxgnygnz *= ex->gnxgny; 978f4bdf6c4SBarry Smith ierr = DMDAGetCorners(ex->da,&ex->xs,&ex->ys,&ex->zs,&ex->nx,&ex->ny,&ex->nz);CHKERRQ(ierr); 979f4bdf6c4SBarry Smith ex->nxny = ex->nx*ex->ny; 980f4bdf6c4SBarry Smith ex->nxnynz = ex->nz*ex->nxny; 981f4bdf6c4SBarry Smith PetscFunctionReturn(0); 982f4bdf6c4SBarry Smith } 983f4bdf6c4SBarry Smith 984f4bdf6c4SBarry Smith #undef __FUNCT__ 985f4bdf6c4SBarry Smith #define __FUNCT__ "MatMult_HYPRESStruct" 986f4bdf6c4SBarry Smith PetscErrorCode MatMult_HYPRESStruct(Mat A,Vec x,Vec y) 987f4bdf6c4SBarry Smith { 988f4bdf6c4SBarry Smith PetscErrorCode ierr; 989f4bdf6c4SBarry Smith PetscScalar *xx,*yy; 990f4bdf6c4SBarry Smith int ilower[3],iupper[3]; 991f4bdf6c4SBarry Smith Mat_HYPRESStruct *mx = (Mat_HYPRESStruct *)(A->data); 992f4bdf6c4SBarry Smith int ordering= mx->dofs_order; 993f4bdf6c4SBarry Smith int nvars= mx->nvars; 994f4bdf6c4SBarry Smith int part= 0; 995f4bdf6c4SBarry Smith int size; 996f4bdf6c4SBarry Smith int i; 997f4bdf6c4SBarry Smith 998f4bdf6c4SBarry Smith PetscFunctionBegin; 999f4bdf6c4SBarry Smith ierr = DMDAGetCorners(mx->da,&ilower[0],&ilower[1],&ilower[2],&iupper[0],&iupper[1],&iupper[2]);CHKERRQ(ierr); 1000f4bdf6c4SBarry Smith iupper[0] += ilower[0] - 1; 1001f4bdf6c4SBarry Smith iupper[1] += ilower[1] - 1; 1002f4bdf6c4SBarry Smith iupper[2] += ilower[2] - 1; 1003f4bdf6c4SBarry Smith 1004f4bdf6c4SBarry Smith size= 1; 1005f4bdf6c4SBarry Smith for (i= 0; i< 3; i++) { 1006f4bdf6c4SBarry Smith size*= (iupper[i]-ilower[i]+1); 1007f4bdf6c4SBarry Smith } 1008f4bdf6c4SBarry Smith 1009f4bdf6c4SBarry Smith /* copy x values over to hypre for variable ordering */ 1010f4bdf6c4SBarry Smith if (ordering) { 1011f4bdf6c4SBarry Smith PetscStackCallHypre(0,HYPRE_SStructVectorSetConstantValues,(mx->ss_b,0.0)); 1012f4bdf6c4SBarry Smith ierr = VecGetArray(x,&xx);CHKERRQ(ierr); 1013f4bdf6c4SBarry Smith for (i= 0; i< nvars; i++) { 1014f4bdf6c4SBarry Smith PetscStackCallHypre(0,HYPRE_SStructVectorSetBoxValues,(mx->ss_b,part,ilower,iupper,i,xx+(size*i))); 1015f4bdf6c4SBarry Smith } 1016f4bdf6c4SBarry Smith ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr); 1017f4bdf6c4SBarry Smith PetscStackCallHypre(0,HYPRE_SStructVectorAssemble,(mx->ss_b)); 1018f4bdf6c4SBarry Smith PetscStackCallHypre(0,HYPRE_SStructMatrixMatvec,(1.0,mx->ss_mat,mx->ss_b,0.0,mx->ss_x)); 1019f4bdf6c4SBarry Smith 1020f4bdf6c4SBarry Smith /* copy solution values back to PETSc */ 1021f4bdf6c4SBarry Smith ierr = VecGetArray(y,&yy);CHKERRQ(ierr); 1022f4bdf6c4SBarry Smith for (i= 0; i< nvars; i++) { 1023f4bdf6c4SBarry Smith PetscStackCallHypre(0,HYPRE_SStructVectorGetBoxValues,(mx->ss_x,part,ilower,iupper,i,yy+(size*i))); 1024f4bdf6c4SBarry Smith } 1025f4bdf6c4SBarry Smith ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr); 1026f4bdf6c4SBarry Smith } else { /* nodal ordering must be mapped to variable ordering for sys_pfmg */ 1027f4bdf6c4SBarry Smith PetscScalar *z; 1028f4bdf6c4SBarry Smith int j, k; 1029f4bdf6c4SBarry Smith 1030f4bdf6c4SBarry Smith ierr = PetscMalloc(nvars*size*sizeof(PetscScalar),&z);CHKERRQ(ierr); 1031f4bdf6c4SBarry Smith PetscStackCallHypre(0,HYPRE_SStructVectorSetConstantValues,(mx->ss_b,0.0)); 1032f4bdf6c4SBarry Smith ierr = VecGetArray(x,&xx);CHKERRQ(ierr); 1033f4bdf6c4SBarry Smith 1034f4bdf6c4SBarry Smith /* transform nodal to hypre's variable ordering for sys_pfmg */ 1035f4bdf6c4SBarry Smith for (i= 0; i< size; i++) { 1036f4bdf6c4SBarry Smith k= i*nvars; 1037f4bdf6c4SBarry Smith for (j= 0; j< nvars; j++) { 1038f4bdf6c4SBarry Smith z[j*size+i]= xx[k+j]; 1039f4bdf6c4SBarry Smith } 1040f4bdf6c4SBarry Smith } 1041f4bdf6c4SBarry Smith for (i= 0; i< nvars; i++) { 1042f4bdf6c4SBarry Smith PetscStackCallHypre(0,HYPRE_SStructVectorSetBoxValues,(mx->ss_b,part,ilower,iupper,i,z+(size*i))); 1043f4bdf6c4SBarry Smith } 1044f4bdf6c4SBarry Smith ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr); 1045f4bdf6c4SBarry Smith PetscStackCallHypre(0,HYPRE_SStructVectorAssemble,(mx->ss_b)); 1046f4bdf6c4SBarry Smith PetscStackCallHypre(0,HYPRE_SStructMatrixMatvec,(1.0,mx->ss_mat,mx->ss_b,0.0,mx->ss_x)); 1047f4bdf6c4SBarry Smith 1048f4bdf6c4SBarry Smith /* copy solution values back to PETSc */ 1049f4bdf6c4SBarry Smith ierr = VecGetArray(y,&yy);CHKERRQ(ierr); 1050f4bdf6c4SBarry Smith for (i= 0; i< nvars; i++) { 1051f4bdf6c4SBarry Smith PetscStackCallHypre(0,HYPRE_SStructVectorGetBoxValues,(mx->ss_x,part,ilower,iupper,i,z+(size*i))); 1052f4bdf6c4SBarry Smith } 1053f4bdf6c4SBarry Smith /* transform hypre's variable ordering for sys_pfmg to nodal ordering */ 1054f4bdf6c4SBarry Smith for (i= 0; i< size; i++) { 1055f4bdf6c4SBarry Smith k= i*nvars; 1056f4bdf6c4SBarry Smith for (j= 0; j< nvars; j++) { 1057f4bdf6c4SBarry Smith yy[k+j]= z[j*size+i]; 1058f4bdf6c4SBarry Smith } 1059f4bdf6c4SBarry Smith } 1060f4bdf6c4SBarry Smith ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr); 1061f4bdf6c4SBarry Smith ierr = PetscFree(z);CHKERRQ(ierr); 1062f4bdf6c4SBarry Smith } 1063f4bdf6c4SBarry Smith PetscFunctionReturn(0); 1064f4bdf6c4SBarry Smith } 1065f4bdf6c4SBarry Smith 1066f4bdf6c4SBarry Smith #undef __FUNCT__ 1067f4bdf6c4SBarry Smith #define __FUNCT__ "MatAssemblyEnd_HYPRESStruct" 1068f4bdf6c4SBarry Smith PetscErrorCode MatAssemblyEnd_HYPRESStruct(Mat mat,MatAssemblyType mode) 1069f4bdf6c4SBarry Smith { 1070f4bdf6c4SBarry Smith Mat_HYPRESStruct *ex = (Mat_HYPRESStruct*) mat->data; 1071f4bdf6c4SBarry Smith PetscErrorCode ierr; 1072f4bdf6c4SBarry Smith 1073f4bdf6c4SBarry Smith PetscFunctionBegin; 1074f4bdf6c4SBarry Smith PetscStackCallHypre(0,HYPRE_SStructMatrixAssemble,(ex->ss_mat)); 1075f4bdf6c4SBarry Smith PetscFunctionReturn(0); 1076f4bdf6c4SBarry Smith } 1077f4bdf6c4SBarry Smith 1078f4bdf6c4SBarry Smith #undef __FUNCT__ 1079f4bdf6c4SBarry Smith #define __FUNCT__ "MatZeroEntries_HYPRESStruct" 1080f4bdf6c4SBarry Smith PetscErrorCode MatZeroEntries_HYPRESStruct(Mat mat) 1081f4bdf6c4SBarry Smith { 1082f4bdf6c4SBarry Smith PetscFunctionBegin; 1083f4bdf6c4SBarry Smith /* before the DMDA is set to the matrix the zero doesn't need to do anything */ 1084f4bdf6c4SBarry Smith PetscFunctionReturn(0); 1085f4bdf6c4SBarry Smith } 1086f4bdf6c4SBarry Smith 1087f4bdf6c4SBarry Smith 1088f4bdf6c4SBarry Smith #undef __FUNCT__ 1089f4bdf6c4SBarry Smith #define __FUNCT__ "MatDestroy_HYPRESStruct" 1090f4bdf6c4SBarry Smith PetscErrorCode MatDestroy_HYPRESStruct(Mat mat) 1091f4bdf6c4SBarry Smith { 1092f4bdf6c4SBarry Smith Mat_HYPRESStruct *ex = (Mat_HYPRESStruct*) mat->data; 1093f4bdf6c4SBarry Smith PetscErrorCode ierr; 1094f4bdf6c4SBarry Smith 1095f4bdf6c4SBarry Smith PetscFunctionBegin; 1096f4bdf6c4SBarry Smith PetscStackCallHypre(0,HYPRE_SStructGraphDestroy,(ex->ss_graph)); 1097f4bdf6c4SBarry Smith PetscStackCallHypre(0,HYPRE_SStructMatrixDestroy,(ex->ss_mat)); 1098f4bdf6c4SBarry Smith PetscStackCallHypre(0,HYPRE_SStructVectorDestroy,(ex->ss_x)); 1099f4bdf6c4SBarry Smith PetscStackCallHypre(0,HYPRE_SStructVectorDestroy,(ex->ss_b)); 1100f4bdf6c4SBarry Smith PetscFunctionReturn(0); 1101f4bdf6c4SBarry Smith } 1102f4bdf6c4SBarry Smith 1103f4bdf6c4SBarry Smith EXTERN_C_BEGIN 1104f4bdf6c4SBarry Smith #undef __FUNCT__ 1105f4bdf6c4SBarry Smith #define __FUNCT__ "MatCreate_HYPRESStruct" 1106f4bdf6c4SBarry Smith PetscErrorCode MatCreate_HYPRESStruct(Mat B) 1107f4bdf6c4SBarry Smith { 1108f4bdf6c4SBarry Smith Mat_HYPRESStruct *ex; 1109f4bdf6c4SBarry Smith PetscErrorCode ierr; 1110f4bdf6c4SBarry Smith 1111f4bdf6c4SBarry Smith PetscFunctionBegin; 1112f4bdf6c4SBarry Smith ierr = PetscNewLog(B,Mat_HYPRESStruct,&ex);CHKERRQ(ierr); 1113f4bdf6c4SBarry Smith B->data = (void*)ex; 1114f4bdf6c4SBarry Smith B->rmap->bs = 1; 1115f4bdf6c4SBarry Smith B->assembled = PETSC_FALSE; 1116f4bdf6c4SBarry Smith 1117f4bdf6c4SBarry Smith B->insertmode = NOT_SET_VALUES; 1118f4bdf6c4SBarry Smith 1119f4bdf6c4SBarry Smith B->ops->assemblyend = MatAssemblyEnd_HYPRESStruct; 1120f4bdf6c4SBarry Smith B->ops->mult = MatMult_HYPRESStruct; 1121f4bdf6c4SBarry Smith B->ops->zeroentries = MatZeroEntries_HYPRESStruct; 1122f4bdf6c4SBarry Smith B->ops->destroy = MatDestroy_HYPRESStruct; 1123f4bdf6c4SBarry Smith 1124f4bdf6c4SBarry Smith ex->needsinitialization = PETSC_TRUE; 1125f4bdf6c4SBarry Smith 1126f4bdf6c4SBarry Smith ierr = MPI_Comm_dup(((PetscObject)B)->comm,&(ex->hcomm));CHKERRQ(ierr); 112795ee5b0eSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSetDM_C","MatSetDM_HYPRESStruct",MatSetDM_HYPRESStruct);CHKERRQ(ierr); 1128f4bdf6c4SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)B,MATHYPRESSTRUCT);CHKERRQ(ierr); 1129f4bdf6c4SBarry Smith PetscFunctionReturn(0); 1130f4bdf6c4SBarry Smith } 1131f4bdf6c4SBarry Smith EXTERN_C_END 1132f4bdf6c4SBarry Smith 1133f4bdf6c4SBarry Smith 1134f4bdf6c4SBarry Smith 1135