Actual source code: mhyp.c
1: #define PETSCMAT_DLL
3: /*
4: Creates hypre ijmatrix from PETSc matrix
5: */
7: #include include/private/matimpl.h
9: #include "HYPRE.h"
10: #include "HYPRE_parcsr_ls.h"
15: PetscErrorCode MatHYPRE_IJMatrixPreallocate(Mat A_d, Mat A_o,HYPRE_IJMatrix ij)
16: {
18: PetscInt i;
19: PetscInt n_d,*ia_d,*ja_d,n_o,*ia_o,*ja_o;
20: PetscTruth done_d=PETSC_FALSE,done_o=PETSC_FALSE;
21: PetscInt *nnz_d=PETSC_NULL,*nnz_o=PETSC_NULL;
22:
24: if (A_d) { /* determine number of nonzero entries in local diagonal part */
25: MatGetRowIJ(A_d,0,PETSC_FALSE,PETSC_FALSE,&n_d,&ia_d,&ja_d,&done_d);
26: if (done_d) {
27: PetscMalloc(n_d*sizeof(PetscInt),&nnz_d);
28: for (i=0; i<n_d; i++) {
29: nnz_d[i] = ia_d[i+1] - ia_d[i];
30: }
31: }
32: MatRestoreRowIJ(A_d,0,PETSC_FALSE,PETSC_FALSE,&n_d,&ia_d,&ja_d,&done_d);
33: }
34: if (A_o) { /* determine number of nonzero entries in local off-diagonal part */
35: MatGetRowIJ(A_o,0,PETSC_FALSE,PETSC_FALSE,&n_o,&ia_o,&ja_o,&done_o);
36: if (done_o) {
37: PetscMalloc(n_o*sizeof(PetscInt),&nnz_o);
38: for (i=0; i<n_o; i++) {
39: nnz_o[i] = ia_o[i+1] - ia_o[i];
40: }
41: }
42: MatRestoreRowIJ(A_o,0,PETSC_FALSE,PETSC_FALSE,&n_o,&ia_o,&ja_o,&done_o);
43: }
44: if (done_d) { /* set number of nonzeros in HYPRE IJ matrix */
45: if (!done_o) { /* only diagonal part */
46: HYPRE_IJMatrixSetRowSizes(ij,nnz_d);
47: } else { /* diagonal and off-diagonal part */
48: HYPRE_IJMatrixSetDiagOffdSizes(ij,nnz_d,nnz_o);
49: }
50: PetscFree(nnz_d);
51: PetscFree(nnz_o);
52: }
53: return(0);
54: }
59: PetscErrorCode MatHYPRE_IJMatrixCreate(Mat A,HYPRE_IJMatrix *ij)
60: {
62: int rstart,rend,cstart,cend;
63:
68: MatPreallocated(A);
69: rstart = A->rmap.rstart;
70: rend = A->rmap.rend;
71: cstart = A->cmap.rstart;
72: cend = A->cmap.rend;
73: HYPRE_IJMatrixCreate(A->comm,rstart,rend-1,cstart,cend-1,ij);
74: HYPRE_IJMatrixSetObjectType(*ij,HYPRE_PARCSR);
75: {
76: PetscTruth same;
77: Mat A_d,A_o;
78: PetscInt *colmap;
79: PetscTypeCompare((PetscObject)A,MATMPIAIJ,&same);
80: if (same) {
81: MatMPIAIJGetSeqAIJ(A,&A_d,&A_o,&colmap);
82: MatHYPRE_IJMatrixPreallocate(A_d,A_o,*ij);
83: return(0);
84: }
85: PetscTypeCompare((PetscObject)A,MATMPIBAIJ,&same);
86: if (same) {
87: MatMPIBAIJGetSeqBAIJ(A,&A_d,&A_o,&colmap);
88: MatHYPRE_IJMatrixPreallocate(A_d,A_o,*ij);
89: return(0);
90: }
91: PetscTypeCompare((PetscObject)A,MATSEQAIJ,&same);
92: if (same) {
93: MatHYPRE_IJMatrixPreallocate(A,PETSC_NULL,*ij);
94: return(0);
95: }
96: PetscTypeCompare((PetscObject)A,MATSEQBAIJ,&same);
97: if (same) {
98: MatHYPRE_IJMatrixPreallocate(A,PETSC_NULL,*ij);
99: return(0);
100: }
101: }
102: return(0);
103: }
105: /*
106: Currently this only works with the MPIAIJ PETSc matrices to make
107: the conversion efficient
108: */
112: PetscErrorCode MatHYPRE_IJMatrixCopy(Mat A,HYPRE_IJMatrix ij)
113: {
114: PetscErrorCode ierr;
115: int i,rstart,rend,ncols;
116: const PetscScalar *values;
117: const int *cols;
121: HYPRE_IJMatrixInitialize(ij);
122: MatGetOwnershipRange(A,&rstart,&rend);
123: for (i=rstart; i<rend; i++) {
124: MatGetRow(A,i,&ncols,&cols,&values);
125: HYPRE_IJMatrixSetValues(ij,1,&ncols,&i,cols,values);
126: MatRestoreRow(A,i,&ncols,&cols,&values);
127: }
128: HYPRE_IJMatrixAssemble(ij);
130: return(0);
131: }