1 2 /* 3 Creates hypre ijmatrix from PETSc matrix 4 */ 5 #include <petscsys.h> 6 #include <petsc/private/matimpl.h> 7 #include <../src/mat/impls/hypre/mhypre.h> 8 #include <../src/mat/impls/aij/mpi/mpiaij.h> 9 #include <HYPRE_struct_mv.h> 10 #include <HYPRE_struct_ls.h> 11 #include <_hypre_struct_mv.h> 12 13 static PetscErrorCode MatHYPRE_CreateFromMat(Mat,Mat_HYPRE*); 14 static PetscErrorCode MatHYPRE_IJMatrixPreallocate(Mat,Mat,HYPRE_IJMatrix); 15 static PetscErrorCode MatHYPRE_IJMatrixFastCopy_MPIAIJ(Mat,HYPRE_IJMatrix); 16 static PetscErrorCode MatHYPRE_IJMatrixFastCopy_SeqAIJ(Mat,HYPRE_IJMatrix); 17 static PetscErrorCode MatHYPRE_MultKernel_Private(Mat,Vec,Vec,PetscBool); 18 19 #undef __FUNCT__ 20 #define __FUNCT__ "MatHYPRE_IJMatrixPreallocate" 21 static PetscErrorCode MatHYPRE_IJMatrixPreallocate(Mat A_d, Mat A_o, HYPRE_IJMatrix ij) 22 { 23 PetscErrorCode ierr; 24 PetscInt i,n_d,n_o; 25 const PetscInt *ia_d,*ia_o; 26 PetscBool done_d=PETSC_FALSE,done_o=PETSC_FALSE; 27 PetscInt *nnz_d=NULL,*nnz_o=NULL; 28 29 PetscFunctionBegin; 30 if (A_d) { /* determine number of nonzero entries in local diagonal part */ 31 ierr = MatGetRowIJ(A_d,0,PETSC_FALSE,PETSC_FALSE,&n_d,&ia_d,NULL,&done_d);CHKERRQ(ierr); 32 if (done_d) { 33 ierr = PetscMalloc1(n_d,&nnz_d);CHKERRQ(ierr); 34 for (i=0; i<n_d; i++) { 35 nnz_d[i] = ia_d[i+1] - ia_d[i]; 36 } 37 } 38 ierr = MatRestoreRowIJ(A_d,0,PETSC_FALSE,PETSC_FALSE,NULL,&ia_d,NULL,&done_d);CHKERRQ(ierr); 39 } 40 if (A_o) { /* determine number of nonzero entries in local off-diagonal part */ 41 ierr = MatGetRowIJ(A_o,0,PETSC_FALSE,PETSC_FALSE,&n_o,&ia_o,NULL,&done_o);CHKERRQ(ierr); 42 if (done_o) { 43 ierr = PetscMalloc1(n_o,&nnz_o);CHKERRQ(ierr); 44 for (i=0; i<n_o; i++) { 45 nnz_o[i] = ia_o[i+1] - ia_o[i]; 46 } 47 } 48 ierr = MatRestoreRowIJ(A_o,0,PETSC_FALSE,PETSC_FALSE,&n_o,&ia_o,NULL,&done_o);CHKERRQ(ierr); 49 } 50 if (done_d) { /* set number of nonzeros in HYPRE IJ matrix */ 51 if (!done_o) { /* only diagonal part */ 52 ierr = PetscMalloc1(n_d,&nnz_o);CHKERRQ(ierr); 53 for (i=0; i<n_d; i++) { 54 nnz_o[i] = 0; 55 } 56 } 57 PetscStackCallStandard(HYPRE_IJMatrixSetDiagOffdSizes,(ij,(HYPRE_Int *)nnz_d,(HYPRE_Int *)nnz_o)); 58 ierr = PetscFree(nnz_d);CHKERRQ(ierr); 59 ierr = PetscFree(nnz_o);CHKERRQ(ierr); 60 } 61 PetscFunctionReturn(0); 62 } 63 64 #undef __FUNCT__ 65 #define __FUNCT__ "MatHYPRE_CreateFromMat" 66 static PetscErrorCode MatHYPRE_CreateFromMat(Mat A, Mat_HYPRE *hA) 67 { 68 PetscErrorCode ierr; 69 PetscInt rstart,rend,cstart,cend; 70 71 PetscFunctionBegin; 72 ierr = MatSetUp(A);CHKERRQ(ierr); 73 rstart = A->rmap->rstart; 74 rend = A->rmap->rend; 75 cstart = A->cmap->rstart; 76 cend = A->cmap->rend; 77 PetscStackCallStandard(HYPRE_IJMatrixCreate,(hA->comm,rstart,rend-1,cstart,cend-1,&hA->ij)); 78 PetscStackCallStandard(HYPRE_IJMatrixSetObjectType,(hA->ij,HYPRE_PARCSR)); 79 { 80 PetscBool same; 81 Mat A_d,A_o; 82 const PetscInt *colmap; 83 ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&same);CHKERRQ(ierr); 84 if (same) { 85 ierr = MatMPIAIJGetSeqAIJ(A,&A_d,&A_o,&colmap);CHKERRQ(ierr); 86 ierr = MatHYPRE_IJMatrixPreallocate(A_d,A_o,hA->ij);CHKERRQ(ierr); 87 PetscFunctionReturn(0); 88 } 89 ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIBAIJ,&same);CHKERRQ(ierr); 90 if (same) { 91 ierr = MatMPIBAIJGetSeqBAIJ(A,&A_d,&A_o,&colmap);CHKERRQ(ierr); 92 ierr = MatHYPRE_IJMatrixPreallocate(A_d,A_o,hA->ij);CHKERRQ(ierr); 93 PetscFunctionReturn(0); 94 } 95 ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&same);CHKERRQ(ierr); 96 if (same) { 97 ierr = MatHYPRE_IJMatrixPreallocate(A,NULL,hA->ij);CHKERRQ(ierr); 98 PetscFunctionReturn(0); 99 } 100 ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQBAIJ,&same);CHKERRQ(ierr); 101 if (same) { 102 ierr = MatHYPRE_IJMatrixPreallocate(A,NULL,hA->ij);CHKERRQ(ierr); 103 PetscFunctionReturn(0); 104 } 105 } 106 PetscFunctionReturn(0); 107 } 108 109 #undef __FUNCT__ 110 #define __FUNCT__ "MatHYPRE_IJMatrixCopy" 111 static 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 PetscStackCallStandard(HYPRE_IJMatrixInitialize,(ij)); 133 ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr); 134 for (i=rstart; i<rend; i++) { 135 ierr = MatGetRow(A,i,&ncols,&cols,&values);CHKERRQ(ierr); 136 PetscStackCallStandard(HYPRE_IJMatrixSetValues,(ij,1,(HYPRE_Int *)&ncols,(HYPRE_Int *)&i,(HYPRE_Int *)cols,values)); 137 ierr = MatRestoreRow(A,i,&ncols,&cols,&values);CHKERRQ(ierr); 138 } 139 PetscFunctionReturn(0); 140 } 141 142 #undef __FUNCT__ 143 #define __FUNCT__ "MatHYPRE_IJMatrixFastCopy_SeqAIJ" 144 static PetscErrorCode MatHYPRE_IJMatrixFastCopy_SeqAIJ(Mat A, HYPRE_IJMatrix ij) 145 { 146 PetscErrorCode ierr; 147 Mat_SeqAIJ *pdiag = (Mat_SeqAIJ*)A->data; 148 149 hypre_ParCSRMatrix *par_matrix; 150 hypre_AuxParCSRMatrix *aux_matrix; 151 hypre_CSRMatrix *hdiag; 152 153 PetscFunctionBegin; 154 PetscStackCallStandard(HYPRE_IJMatrixInitialize,(ij)); 155 par_matrix = (hypre_ParCSRMatrix*)hypre_IJMatrixObject(ij); 156 aux_matrix = (hypre_AuxParCSRMatrix*)hypre_IJMatrixTranslator(ij); 157 hdiag = hypre_ParCSRMatrixDiag(par_matrix); 158 /* 159 this is the Hack part where we monkey directly with the hypre datastructures 160 */ 161 ierr = PetscMemcpy(hdiag->i,pdiag->i,(A->rmap->n + 1)*sizeof(PetscInt)); 162 ierr = PetscMemcpy(hdiag->j,pdiag->j,pdiag->nz*sizeof(PetscInt)); 163 ierr = PetscMemcpy(hdiag->data,pdiag->a,pdiag->nz*sizeof(PetscScalar)); 164 hypre_AuxParCSRMatrixNeedAux(aux_matrix) = 0; 165 PetscFunctionReturn(0); 166 } 167 168 #undef __FUNCT__ 169 #define __FUNCT__ "MatHYPRE_IJMatrixFastCopy_MPIAIJ" 170 static PetscErrorCode MatHYPRE_IJMatrixFastCopy_MPIAIJ(Mat A, HYPRE_IJMatrix ij) 171 { 172 PetscErrorCode ierr; 173 Mat_MPIAIJ *pA = (Mat_MPIAIJ*)A->data; 174 Mat_SeqAIJ *pdiag,*poffd; 175 PetscInt i,*garray = pA->garray,*jj,cstart,*pjj; 176 177 hypre_ParCSRMatrix *par_matrix; 178 hypre_AuxParCSRMatrix *aux_matrix; 179 hypre_CSRMatrix *hdiag,*hoffd; 180 181 PetscFunctionBegin; 182 pdiag = (Mat_SeqAIJ*) pA->A->data; 183 poffd = (Mat_SeqAIJ*) pA->B->data; 184 /* cstart is only valid for square MPIAIJ layed out in the usual way */ 185 ierr = MatGetOwnershipRange(A,&cstart,NULL);CHKERRQ(ierr); 186 187 PetscStackCallStandard(HYPRE_IJMatrixInitialize,(ij)); 188 par_matrix = (hypre_ParCSRMatrix*)hypre_IJMatrixObject(ij); 189 aux_matrix = (hypre_AuxParCSRMatrix*)hypre_IJMatrixTranslator(ij); 190 hdiag = hypre_ParCSRMatrixDiag(par_matrix); 191 hoffd = hypre_ParCSRMatrixOffd(par_matrix); 192 193 /* 194 this is the Hack part where we monkey directly with the hypre datastructures 195 */ 196 ierr = PetscMemcpy(hdiag->i,pdiag->i,(pA->A->rmap->n + 1)*sizeof(PetscInt)); 197 /* need to shift the diag column indices (hdiag->j) back to global numbering since hypre is expecting this */ 198 jj = (PetscInt*)hdiag->j; 199 pjj = (PetscInt*)pdiag->j; 200 for (i=0; i<pdiag->nz; i++) jj[i] = cstart + pjj[i]; 201 ierr = PetscMemcpy(hdiag->data,pdiag->a,pdiag->nz*sizeof(PetscScalar)); 202 ierr = PetscMemcpy(hoffd->i,poffd->i,(pA->A->rmap->n + 1)*sizeof(PetscInt)); 203 /* need to move the offd column indices (hoffd->j) back to global numbering since hypre is expecting this 204 If we hacked a hypre a bit more we might be able to avoid this step */ 205 jj = (PetscInt*) hoffd->j; 206 pjj = (PetscInt*) poffd->j; 207 for (i=0; i<poffd->nz; i++) jj[i] = garray[pjj[i]]; 208 ierr = PetscMemcpy(hoffd->data,poffd->a,poffd->nz*sizeof(PetscScalar)); 209 210 hypre_AuxParCSRMatrixNeedAux(aux_matrix) = 0; 211 PetscFunctionReturn(0); 212 } 213 214 #undef __FUNCT__ 215 #define __FUNCT__ "MatConvert_AIJ_HYPRE" 216 PETSC_INTERN PetscErrorCode MatConvert_AIJ_HYPRE(Mat A, MatType type, MatReuse reuse, Mat *B) 217 { 218 Mat_HYPRE *hB; 219 MPI_Comm comm = PetscObjectComm((PetscObject)A); 220 PetscErrorCode ierr; 221 222 PetscFunctionBegin; 223 if (reuse == MAT_INPLACE_MATRIX) SETERRQ(comm,PETSC_ERR_SUP,"Unsupported MAT_INPLACE_MATRIX"); 224 if (reuse == MAT_REUSE_MATRIX) { 225 /* always destroy the old matrix and create a new memory; 226 hope this does not churn the memory too much. The problem 227 is I do not know if it is possible to put the matrix back to 228 its initial state so that we can directly copy the values 229 the second time through. */ 230 hB = (Mat_HYPRE*)((*B)->data); 231 PetscStackCallStandard(HYPRE_IJMatrixDestroy,(hB->ij)); 232 } else { 233 ierr = MatCreate(comm,B);CHKERRQ(ierr); 234 ierr = MatSetType(*B,MATHYPRE);CHKERRQ(ierr); 235 ierr = MatSetSizes(*B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr); 236 ierr = MatSetUp(*B);CHKERRQ(ierr); 237 hB = (Mat_HYPRE*)((*B)->data); 238 } 239 ierr = MatHYPRE_CreateFromMat(A,hB);CHKERRQ(ierr); 240 ierr = MatHYPRE_IJMatrixCopy(A,hB->ij);CHKERRQ(ierr); 241 ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 242 ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 243 PetscFunctionReturn(0); 244 } 245 246 #undef __FUNCT__ 247 #define __FUNCT__ "MatConvert_HYPRE_AIJ" 248 PetscErrorCode MatConvert_HYPRE_AIJ(Mat A, MatType mtype, MatReuse reuse, Mat *B) 249 { 250 Mat_HYPRE *hA = (Mat_HYPRE*)A->data; 251 hypre_ParCSRMatrix *parcsr; 252 hypre_CSRMatrix *hdiag,*hoffd; 253 MPI_Comm comm; 254 PetscScalar *da,*oa,*aptr; 255 PetscInt *dii,*djj,*oii,*ojj,*iptr; 256 PetscInt i,dnnz,onnz,m,n; 257 int type; 258 PetscMPIInt size; 259 PetscErrorCode ierr; 260 261 PetscFunctionBegin; 262 comm = PetscObjectComm((PetscObject)A); 263 PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(hA->ij,&type)); 264 if (type != HYPRE_PARCSR) SETERRQ(comm,PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported"); 265 if (reuse == MAT_REUSE_MATRIX) { 266 PetscBool ismpiaij,isseqaij; 267 ierr = PetscObjectTypeCompare((PetscObject)*B,MATMPIAIJ,&ismpiaij);CHKERRQ(ierr); 268 ierr = PetscObjectTypeCompare((PetscObject)*B,MATSEQAIJ,&isseqaij);CHKERRQ(ierr); 269 if (!ismpiaij && !isseqaij) SETERRQ(comm,PETSC_ERR_SUP,"Only MATMPIAIJ or MATSEQAIJ are supported"); 270 } 271 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 272 273 PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hA->ij,(void**)&parcsr)); 274 hdiag = hypre_ParCSRMatrixDiag(parcsr); 275 hoffd = hypre_ParCSRMatrixOffd(parcsr); 276 m = hypre_CSRMatrixNumRows(hdiag); 277 n = hypre_CSRMatrixNumCols(hdiag); 278 dnnz = hypre_CSRMatrixNumNonzeros(hdiag); 279 onnz = hypre_CSRMatrixNumNonzeros(hoffd); 280 if (reuse != MAT_REUSE_MATRIX) { 281 ierr = PetscMalloc1(m+1,&dii);CHKERRQ(ierr); 282 ierr = PetscMalloc1(dnnz,&djj);CHKERRQ(ierr); 283 ierr = PetscMalloc1(dnnz,&da);CHKERRQ(ierr); 284 } else { 285 PetscInt nr; 286 PetscBool done; 287 if (size > 1) { 288 Mat_MPIAIJ *b = (Mat_MPIAIJ*)((*B)->data); 289 290 ierr = MatGetRowIJ(b->A,0,PETSC_FALSE,PETSC_FALSE,&nr,(const PetscInt**)&dii,(const PetscInt**)&djj,&done);CHKERRQ(ierr); 291 if (nr != m) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot reuse mat: invalid number of local rows in diag part! %D != %D",nr,m); 292 if (dii[nr] < dnnz) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot reuse mat: invalid number of nonzeros in diag part! reuse %D hypre %D",dii[nr],dnnz); 293 ierr = MatSeqAIJGetArray(b->A,&da);CHKERRQ(ierr); 294 } else { 295 ierr = MatGetRowIJ(*B,0,PETSC_FALSE,PETSC_FALSE,&nr,(const PetscInt**)&dii,(const PetscInt**)&djj,&done);CHKERRQ(ierr); 296 if (nr != m) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot reuse mat: invalid number of local rows! %D != %D",nr,m); 297 if (dii[nr] < dnnz) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot reuse mat: invalid number of nonzeros! reuse %D hypre %D",dii[nr],dnnz); 298 ierr = MatSeqAIJGetArray(*B,&da);CHKERRQ(ierr); 299 } 300 } 301 ierr = PetscMemcpy(dii,hypre_CSRMatrixI(hdiag),(m+1)*sizeof(PetscInt));CHKERRQ(ierr); 302 ierr = PetscMemcpy(djj,hypre_CSRMatrixJ(hdiag),dnnz*sizeof(PetscInt));CHKERRQ(ierr); 303 ierr = PetscMemcpy(da,hypre_CSRMatrixData(hdiag),dnnz*sizeof(PetscScalar));CHKERRQ(ierr); 304 iptr = djj; 305 aptr = da; 306 for (i=0; i<m; i++) { 307 PetscInt nc = dii[i+1]-dii[i]; 308 ierr = PetscSortIntWithScalarArray(nc,iptr,aptr);CHKERRQ(ierr); 309 iptr += nc; 310 aptr += nc; 311 } 312 if (size > 1) { 313 PetscInt *offdj,*coffd; 314 315 if (reuse != MAT_REUSE_MATRIX) { 316 ierr = PetscMalloc1(m+1,&oii);CHKERRQ(ierr); 317 ierr = PetscMalloc1(onnz,&ojj);CHKERRQ(ierr); 318 ierr = PetscMalloc1(onnz,&oa);CHKERRQ(ierr); 319 } else { 320 Mat_MPIAIJ *b = (Mat_MPIAIJ*)((*B)->data); 321 PetscInt nr,hr = hypre_CSRMatrixNumRows(hoffd); 322 PetscBool done; 323 324 ierr = MatGetRowIJ(b->B,0,PETSC_FALSE,PETSC_FALSE,&nr,(const PetscInt**)&oii,(const PetscInt**)&ojj,&done);CHKERRQ(ierr); 325 if (nr != hr) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot reuse mat: invalid number of local rows in offdiag part! %D != %D",nr,hr); 326 if (oii[nr] < onnz) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot reuse mat: invalid number of nonzeros in offdiag part! reuse %D hypre %D",oii[nr],onnz); 327 ierr = MatSeqAIJGetArray(b->B,&oa);CHKERRQ(ierr); 328 } 329 ierr = PetscMemcpy(oii,hypre_CSRMatrixI(hoffd),(m+1)*sizeof(PetscInt));CHKERRQ(ierr); 330 offdj = hypre_CSRMatrixJ(hoffd); 331 coffd = hypre_ParCSRMatrixColMapOffd(parcsr); 332 for (i=0; i<onnz; i++) ojj[i] = coffd[offdj[i]]; 333 ierr = PetscMemcpy(oa,hypre_CSRMatrixData(hoffd),onnz*sizeof(PetscScalar));CHKERRQ(ierr); 334 iptr = ojj; 335 aptr = oa; 336 for (i=0; i<m; i++) { 337 PetscInt nc = oii[i+1]-oii[i]; 338 ierr = PetscSortIntWithScalarArray(nc,iptr,aptr);CHKERRQ(ierr); 339 iptr += nc; 340 aptr += nc; 341 } 342 if (reuse != MAT_REUSE_MATRIX) { 343 Mat_MPIAIJ *b; 344 Mat_SeqAIJ *d,*o; 345 ierr = MatCreateMPIAIJWithSplitArrays(comm,m,n,PETSC_DECIDE,PETSC_DECIDE,dii, 346 djj,da,oii,ojj,oa,B);CHKERRQ(ierr); 347 /* hack MPIAIJ */ 348 b = (Mat_MPIAIJ*)((*B)->data); 349 d = (Mat_SeqAIJ*)b->A->data; 350 o = (Mat_SeqAIJ*)b->B->data; 351 d->free_a = PETSC_TRUE; 352 d->free_ij = PETSC_TRUE; 353 o->free_a = PETSC_TRUE; 354 o->free_ij = PETSC_TRUE; 355 } 356 } else if (reuse != MAT_REUSE_MATRIX) { 357 Mat_SeqAIJ* b; 358 ierr = MatCreateSeqAIJWithArrays(comm,m,n,dii,djj,da,B);CHKERRQ(ierr); 359 /* hack SeqAIJ */ 360 b = (Mat_SeqAIJ*)((*B)->data); 361 b->free_a = PETSC_TRUE; 362 b->free_ij = PETSC_TRUE; 363 } 364 if (reuse == MAT_INPLACE_MATRIX) { 365 ierr = MatHeaderReplace(A,B);CHKERRQ(ierr); 366 } 367 PetscFunctionReturn(0); 368 } 369 370 #undef __FUNCT__ 371 #define __FUNCT__ "MatMultTranspose_HYPRE" 372 PetscErrorCode MatMultTranspose_HYPRE(Mat A, Vec x, Vec y) 373 { 374 PetscErrorCode ierr; 375 376 PetscFunctionBegin; 377 ierr = MatHYPRE_MultKernel_Private(A,x,y,PETSC_TRUE);CHKERRQ(ierr); 378 PetscFunctionReturn(0); 379 } 380 381 #undef __FUNCT__ 382 #define __FUNCT__ "MatMult_HYPRE" 383 PetscErrorCode MatMult_HYPRE(Mat A, Vec x, Vec y) 384 { 385 PetscErrorCode ierr; 386 387 PetscFunctionBegin; 388 ierr = MatHYPRE_MultKernel_Private(A,x,y,PETSC_FALSE);CHKERRQ(ierr); 389 PetscFunctionReturn(0); 390 } 391 392 #undef __FUNCT__ 393 #define __FUNCT__ "MatHYPRE_MultKernel_Private" 394 PetscErrorCode MatHYPRE_MultKernel_Private(Mat A, Vec x, Vec y, PetscBool trans) 395 { 396 Mat_HYPRE *hA = (Mat_HYPRE*)A->data; 397 hypre_ParCSRMatrix *parcsr; 398 hypre_ParVector *hx,*hy; 399 PetscScalar *ax,*ay,*sax,*say; 400 PetscErrorCode ierr; 401 402 PetscFunctionBegin; 403 PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hA->ij,(void**)&parcsr)); 404 PetscStackCallStandard(HYPRE_IJVectorGetObject,(hA->x,(void**)&hx)); 405 PetscStackCallStandard(HYPRE_IJVectorGetObject,(hA->b,(void**)&hy)); 406 ierr = VecGetArrayRead(x,(const PetscScalar**)&ax);CHKERRQ(ierr); 407 ierr = VecGetArray(y,&ay);CHKERRQ(ierr); 408 if (trans) { 409 HYPREReplacePointer(hA->x,ay,say); 410 HYPREReplacePointer(hA->b,ax,sax); 411 hypre_ParCSRMatrixMatvecT(1.,parcsr,hy,0.,hx); 412 HYPREReplacePointer(hA->x,say,ay); 413 HYPREReplacePointer(hA->b,sax,ax); 414 } else { 415 HYPREReplacePointer(hA->x,ax,sax); 416 HYPREReplacePointer(hA->b,ay,say); 417 hypre_ParCSRMatrixMatvec(1.,parcsr,hx,0.,hy); 418 HYPREReplacePointer(hA->x,sax,ax); 419 HYPREReplacePointer(hA->b,say,ay); 420 } 421 ierr = VecRestoreArrayRead(x,(const PetscScalar**)&ax);CHKERRQ(ierr); 422 ierr = VecRestoreArray(y,&ay);CHKERRQ(ierr); 423 PetscFunctionReturn(0); 424 } 425 426 #undef __FUNCT__ 427 #define __FUNCT__ "MatDestroy_HYPRE" 428 PetscErrorCode MatDestroy_HYPRE(Mat A) 429 { 430 Mat_HYPRE *hA = (Mat_HYPRE*)A->data; 431 PetscErrorCode ierr; 432 433 PetscFunctionBegin; 434 if (hA->x) PetscStackCallStandard(HYPRE_IJVectorDestroy,(hA->x)); 435 if (hA->b) PetscStackCallStandard(HYPRE_IJVectorDestroy,(hA->b)); 436 if (hA->ij) PetscStackCallStandard(HYPRE_IJMatrixDestroy,(hA->ij)); 437 if (hA->comm) { ierr = MPI_Comm_free(&hA->comm);CHKERRQ(ierr);} 438 ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_hypre_aij_C",NULL);CHKERRQ(ierr); 439 ierr = PetscFree(A->data);CHKERRQ(ierr); 440 PetscFunctionReturn(0); 441 } 442 443 #undef __FUNCT__ 444 #define __FUNCT__ "MatSetUp_HYPRE" 445 PetscErrorCode MatSetUp_HYPRE(Mat A) 446 { 447 Mat_HYPRE *hA = (Mat_HYPRE*)A->data; 448 Vec x,b; 449 PetscErrorCode ierr; 450 451 PetscFunctionBegin; 452 if (hA->x) PetscFunctionReturn(0); 453 ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr); 454 ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr); 455 ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)A),1,A->cmap->n,A->cmap->N,NULL,&x);CHKERRQ(ierr); 456 ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)A),1,A->rmap->n,A->rmap->N,NULL,&b);CHKERRQ(ierr); 457 ierr = VecHYPRE_IJVectorCreate(x,&hA->x);CHKERRQ(ierr); 458 ierr = VecHYPRE_IJVectorCreate(b,&hA->b);CHKERRQ(ierr); 459 ierr = VecDestroy(&x);CHKERRQ(ierr); 460 ierr = VecDestroy(&b);CHKERRQ(ierr); 461 PetscFunctionReturn(0); 462 } 463 464 #undef __FUNCT__ 465 #define __FUNCT__ "MatAssemblyEnd_HYPRE" 466 PetscErrorCode MatAssemblyEnd_HYPRE(Mat A, MatAssemblyType mode) 467 { 468 Mat_HYPRE *hA = (Mat_HYPRE*)A->data; 469 PetscErrorCode ierr; 470 471 PetscFunctionBegin; 472 PetscStackCallStandard(HYPRE_IJMatrixAssemble,(hA->ij)); 473 PetscFunctionReturn(0); 474 } 475 476 /* Is this needed? 477 int type; 478 hypre_ParCSRMatrix *parcsr; 479 PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(hA->ij,&type)); 480 if (type != HYPRE_PARCSR) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Only HYPRE_PARCSR supported"); 481 PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hA->ij,(void**)&parcsr)); 482 PetscStackCallStandard(hypre_MatvecCommPkgCreate,(parcsr)); 483 */ 484 #undef __FUNCT__ 485 #define __FUNCT__ "MatCreate_HYPRE" 486 PETSC_EXTERN PetscErrorCode MatCreate_HYPRE(Mat B) 487 { 488 Mat_HYPRE *hB; 489 PetscErrorCode ierr; 490 491 PetscFunctionBegin; 492 ierr = PetscNewLog(B,&hB);CHKERRQ(ierr); 493 B->data = (void*)hB; 494 B->rmap->bs = 1; 495 B->assembled = PETSC_FALSE; 496 497 B->ops->mult = MatMult_HYPRE; 498 B->ops->multtranspose = MatMultTranspose_HYPRE; 499 B->ops->setup = MatSetUp_HYPRE; 500 B->ops->destroy = MatDestroy_HYPRE; 501 B->ops->assemblyend = MatAssemblyEnd_HYPRE; 502 503 ierr = MPI_Comm_dup(PetscObjectComm((PetscObject)B),&hB->comm);CHKERRQ(ierr); 504 ierr = PetscObjectChangeTypeName((PetscObject)B,MATHYPRE);CHKERRQ(ierr); 505 ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_hypre_aij_C",MatConvert_HYPRE_AIJ);CHKERRQ(ierr); 506 PetscFunctionReturn(0); 507 } 508 509 #if 0 510 /* 511 Does NOT copy the data over, instead uses DIRECTLY the pointers from the PETSc MPIAIJ format 512 513 This is UNFINISHED and does NOT work! The problem is that hypre puts the diagonal entry first 514 which will corrupt the PETSc data structure if we did this. Need a work around to this problem. 515 */ 516 #include <_hypre_IJ_mv.h> 517 #include <HYPRE_IJ_mv.h> 518 519 #undef __FUNCT__ 520 #define __FUNCT__ "MatHYPRE_IJMatrixLink" 521 PetscErrorCode MatHYPRE_IJMatrixLink(Mat A, HYPRE_IJMatrix *ij) 522 { 523 PetscErrorCode ierr; 524 PetscInt rstart,rend,cstart,cend; 525 PetscBool flg; 526 hypre_AuxParCSRMatrix *aux_matrix; 527 528 PetscFunctionBegin; 529 PetscValidHeaderSpecific(A,MAT_CLASSID,1); 530 PetscValidType(A,1); 531 PetscValidPointer(ij,2); 532 ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&flg);CHKERRQ(ierr); 533 if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Can only use with PETSc MPIAIJ matrices"); 534 ierr = MatSetUp(A);CHKERRQ(ierr); 535 536 ierr = PetscLogEventBegin(MAT_Convert,A,0,0,0);CHKERRQ(ierr); 537 rstart = A->rmap->rstart; 538 rend = A->rmap->rend; 539 cstart = A->cmap->rstart; 540 cend = A->cmap->rend; 541 PetscStackCallStandard(HYPRE_IJMatrixCreate,(PetscObjectComm((PetscObject)A),rstart,rend-1,cstart,cend-1,ij)); 542 PetscStackCallStandard(HYPRE_IJMatrixSetObjectType,(*ij,HYPRE_PARCSR)); 543 544 PetscStackCallStandard(HYPRE_IJMatrixInitialize,(*ij)); 545 PetscStackCall("hypre_IJMatrixTranslator",aux_matrix = (hypre_AuxParCSRMatrix*)hypre_IJMatrixTranslator(*ij)); 546 547 hypre_AuxParCSRMatrixNeedAux(aux_matrix) = 0; 548 549 /* this is the Hack part where we monkey directly with the hypre datastructures */ 550 551 PetscStackCallStandard(HYPRE_IJMatrixAssemble,(*ij)); 552 ierr = PetscLogEventEnd(MAT_Convert,A,0,0,0);CHKERRQ(ierr); 553 PetscFunctionReturn(0); 554 } 555 #endif 556