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