hypre_ParCSRMatrix * hypre_ParCSRBlockMatrixCompress( hypre_ParCSRBlockMatrix *matrix ) { MPI_Comm comm = hypre_ParCSRBlockMatrixComm(matrix); hypre_CSRBlockMatrix *diag = hypre_ParCSRBlockMatrixDiag(matrix); hypre_CSRBlockMatrix *offd = hypre_ParCSRBlockMatrixOffd(matrix); HYPRE_Int global_num_rows = hypre_ParCSRBlockMatrixGlobalNumRows(matrix); HYPRE_Int global_num_cols = hypre_ParCSRBlockMatrixGlobalNumCols(matrix); HYPRE_Int *row_starts = hypre_ParCSRBlockMatrixRowStarts(matrix); HYPRE_Int *col_starts = hypre_ParCSRBlockMatrixColStarts(matrix); HYPRE_Int num_cols_offd = hypre_CSRBlockMatrixNumCols(offd); HYPRE_Int num_nonzeros_diag = hypre_CSRBlockMatrixNumNonzeros(diag); HYPRE_Int num_nonzeros_offd = hypre_CSRBlockMatrixNumNonzeros(offd); hypre_ParCSRMatrix *matrix_C; HYPRE_Int i; matrix_C = hypre_ParCSRMatrixCreate(comm, global_num_rows, global_num_cols, row_starts,col_starts,num_cols_offd,num_nonzeros_diag,num_nonzeros_offd); hypre_ParCSRMatrixInitialize(matrix_C); hypre_CSRMatrixDestroy(hypre_ParCSRMatrixDiag(matrix_C)); hypre_ParCSRMatrixDiag(matrix_C) = hypre_CSRBlockMatrixCompress(diag); hypre_CSRMatrixDestroy(hypre_ParCSRMatrixOffd(matrix_C)); hypre_ParCSRMatrixOffd(matrix_C) = hypre_CSRBlockMatrixCompress(offd); for(i = 0; i < num_cols_offd; i++) hypre_ParCSRMatrixColMapOffd(matrix_C)[i] = hypre_ParCSRBlockMatrixColMapOffd(matrix)[i]; return matrix_C; }
HYPRE_Int hypre_ParCSRBlockMatrixDestroy( hypre_ParCSRBlockMatrix *matrix ) { if (matrix) { if ( hypre_ParCSRBlockMatrixOwnsData(matrix) ) { hypre_CSRBlockMatrixDestroy(hypre_ParCSRBlockMatrixDiag(matrix)); hypre_CSRBlockMatrixDestroy(hypre_ParCSRBlockMatrixOffd(matrix)); if (hypre_ParCSRBlockMatrixColMapOffd(matrix)) hypre_TFree(hypre_ParCSRBlockMatrixColMapOffd(matrix)); if (hypre_ParCSRBlockMatrixCommPkg(matrix)) hypre_MatvecCommPkgDestroy(hypre_ParCSRBlockMatrixCommPkg(matrix)); if (hypre_ParCSRBlockMatrixCommPkgT(matrix)) hypre_MatvecCommPkgDestroy(hypre_ParCSRBlockMatrixCommPkgT(matrix)); } if ( hypre_ParCSRBlockMatrixOwnsRowStarts(matrix) ) hypre_TFree(hypre_ParCSRBlockMatrixRowStarts(matrix)); if ( hypre_ParCSRBlockMatrixOwnsColStarts(matrix) ) hypre_TFree(hypre_ParCSRBlockMatrixColStarts(matrix)); if (hypre_ParCSRBlockMatrixAssumedPartition(matrix)) hypre_ParCSRBlockMatrixDestroyAssumedPartition(matrix); hypre_TFree(matrix); } return hypre_error_flag; }
hypre_ParCSRMatrix * hypre_ParCSRBlockMatrixConvertToParCSRMatrix(hypre_ParCSRBlockMatrix *matrix) { MPI_Comm comm = hypre_ParCSRBlockMatrixComm(matrix); hypre_CSRBlockMatrix *diag = hypre_ParCSRBlockMatrixDiag(matrix); hypre_CSRBlockMatrix *offd = hypre_ParCSRBlockMatrixOffd(matrix); HYPRE_Int block_size = hypre_ParCSRBlockMatrixBlockSize(matrix); HYPRE_Int global_num_rows = hypre_ParCSRBlockMatrixGlobalNumRows(matrix); HYPRE_Int global_num_cols = hypre_ParCSRBlockMatrixGlobalNumCols(matrix); HYPRE_Int *row_starts = hypre_ParCSRBlockMatrixRowStarts(matrix); HYPRE_Int *col_starts = hypre_ParCSRBlockMatrixColStarts(matrix); HYPRE_Int num_cols_offd = hypre_CSRBlockMatrixNumCols(offd); HYPRE_Int num_nonzeros_diag = hypre_CSRBlockMatrixNumNonzeros(diag); HYPRE_Int num_nonzeros_offd = hypre_CSRBlockMatrixNumNonzeros(offd); hypre_ParCSRMatrix *matrix_C; HYPRE_Int *matrix_C_row_starts; HYPRE_Int *matrix_C_col_starts; HYPRE_Int *counter, *new_j_map; HYPRE_Int size_j, size_map, index, new_num_cols, removed = 0; HYPRE_Int *offd_j, *col_map_offd, *new_col_map_offd; HYPRE_Int num_procs, i, j; hypre_CSRMatrix *diag_nozeros, *offd_nozeros; hypre_MPI_Comm_size(comm,&num_procs); #ifdef HYPRE_NO_GLOBAL_PARTITION matrix_C_row_starts = hypre_CTAlloc(HYPRE_Int, 2); matrix_C_col_starts = hypre_CTAlloc(HYPRE_Int, 2); for(i = 0; i < 2; i++) { matrix_C_row_starts[i] = row_starts[i]*block_size; matrix_C_col_starts[i] = col_starts[i]*block_size; } #else matrix_C_row_starts = hypre_CTAlloc(HYPRE_Int, num_procs + 1); matrix_C_col_starts = hypre_CTAlloc(HYPRE_Int, num_procs + 1); for(i = 0; i < num_procs + 1; i++) { matrix_C_row_starts[i] = row_starts[i]*block_size; matrix_C_col_starts[i] = col_starts[i]*block_size; } #endif matrix_C = hypre_ParCSRMatrixCreate(comm, global_num_rows*block_size, global_num_cols*block_size, matrix_C_row_starts, matrix_C_col_starts, num_cols_offd*block_size, num_nonzeros_diag*block_size*block_size, num_nonzeros_offd*block_size*block_size); hypre_ParCSRMatrixInitialize(matrix_C); /* DIAG */ hypre_CSRMatrixDestroy(hypre_ParCSRMatrixDiag(matrix_C)); hypre_ParCSRMatrixDiag(matrix_C) = hypre_CSRBlockMatrixConvertToCSRMatrix(diag); /* AB - added to delete zeros */ diag_nozeros = hypre_CSRMatrixDeleteZeros( hypre_ParCSRMatrixDiag(matrix_C), 1e-14); if(diag_nozeros) { hypre_CSRMatrixDestroy(hypre_ParCSRMatrixDiag(matrix_C)); hypre_ParCSRMatrixDiag(matrix_C) = diag_nozeros; } /* OFF-DIAG */ hypre_CSRMatrixDestroy(hypre_ParCSRMatrixOffd(matrix_C)); hypre_ParCSRMatrixOffd(matrix_C) = hypre_CSRBlockMatrixConvertToCSRMatrix(offd); /* AB - added to delete zeros - this just deletes from data and j arrays */ offd_nozeros = hypre_CSRMatrixDeleteZeros( hypre_ParCSRMatrixOffd(matrix_C), 1e-14); if(offd_nozeros) { hypre_CSRMatrixDestroy(hypre_ParCSRMatrixOffd(matrix_C)); hypre_ParCSRMatrixOffd(matrix_C) = offd_nozeros; removed = 1; } /* now convert the col_map_offd */ for (i = 0; i < num_cols_offd; i++) for (j = 0; j < block_size; j++) hypre_ParCSRMatrixColMapOffd(matrix_C)[i*block_size + j] = hypre_ParCSRBlockMatrixColMapOffd(matrix)[i]*block_size + j; /* if we deleted zeros, then it is possible that col_map_offd can be compressed as well - this requires some amount of work that could be skipped... */ if (removed) { size_map = num_cols_offd*block_size; counter = hypre_CTAlloc(HYPRE_Int, size_map); new_j_map = hypre_CTAlloc(HYPRE_Int, size_map); offd_j = hypre_CSRMatrixJ(hypre_ParCSRMatrixOffd(matrix_C)); col_map_offd = hypre_ParCSRMatrixColMapOffd(matrix_C); size_j = hypre_CSRMatrixNumNonzeros(hypre_ParCSRMatrixOffd(matrix_C)); /* mark which off_d entries are found in j */ for (i=0; i < size_j; i++) { counter[offd_j[i]] = 1; } /*now find new numbering for columns (we will delete the cols where counter = 0*/ index = 0; for (i=0; i < size_map; i++) { if (counter[i]) new_j_map[i] = index++; } new_num_cols = index; /* if there are some col entries to remove: */ if (!(index == size_map)) { /* go thru j and adjust entries */ for (i=0; i < size_j; i++) { offd_j[i] = new_j_map[offd_j[i]]; } /*now go thru col map and get rid of non-needed entries */ new_col_map_offd = hypre_CTAlloc(HYPRE_Int, new_num_cols); index = 0; for (i=0; i < size_map; i++) { if (counter[i]) { new_col_map_offd[index++] = col_map_offd[i]; } } /* set the new col map */ hypre_TFree(col_map_offd); hypre_ParCSRMatrixColMapOffd(matrix_C) = new_col_map_offd; /* modify the number of cols */ hypre_CSRMatrixNumCols(hypre_ParCSRMatrixOffd(matrix_C)) = new_num_cols; } hypre_TFree(new_j_map); hypre_TFree(counter); } hypre_ParCSRMatrixSetNumNonzeros( matrix_C ); hypre_ParCSRMatrixSetDNumNonzeros( matrix_C ); /* we will not copy the comm package */ hypre_ParCSRMatrixCommPkg(matrix_C) = NULL; return matrix_C; }
hypre_ParCSRBlockMatrix * hypre_ParCSRBlockMatrixCreate( MPI_Comm comm, HYPRE_Int block_size, HYPRE_Int global_num_rows, HYPRE_Int global_num_cols, HYPRE_Int *row_starts, HYPRE_Int *col_starts, HYPRE_Int num_cols_offd, HYPRE_Int num_nonzeros_diag, HYPRE_Int num_nonzeros_offd ) { hypre_ParCSRBlockMatrix *matrix; HYPRE_Int num_procs, my_id; HYPRE_Int local_num_rows, local_num_cols; HYPRE_Int first_row_index, first_col_diag; matrix = hypre_CTAlloc(hypre_ParCSRBlockMatrix, 1); hypre_MPI_Comm_rank(comm,&my_id); hypre_MPI_Comm_size(comm,&num_procs); if (!row_starts) { #ifdef HYPRE_NO_GLOBAL_PARTITION hypre_GenerateLocalPartitioning(global_num_rows, num_procs, my_id, &row_starts); #else hypre_GeneratePartitioning(global_num_rows,num_procs,&row_starts); #endif } if (!col_starts) { if (global_num_rows == global_num_cols) { col_starts = row_starts; } else { #ifdef HYPRE_NO_GLOBAL_PARTITION hypre_GenerateLocalPartitioning(global_num_cols, num_procs, my_id, &col_starts); #else hypre_GeneratePartitioning(global_num_cols,num_procs,&col_starts); #endif } } #ifdef HYPRE_NO_GLOBAL_PARTITION /* row_starts[0] is start of local rows. row_starts[1] is start of next processor's rows */ first_row_index = row_starts[0]; local_num_rows = row_starts[1]-first_row_index ; first_col_diag = col_starts[0]; local_num_cols = col_starts[1]-first_col_diag; #else first_row_index = row_starts[my_id]; local_num_rows = row_starts[my_id+1]-first_row_index; first_col_diag = col_starts[my_id]; local_num_cols = col_starts[my_id+1]-first_col_diag; #endif hypre_ParCSRBlockMatrixComm(matrix) = comm; hypre_ParCSRBlockMatrixDiag(matrix) = hypre_CSRBlockMatrixCreate(block_size, local_num_rows,local_num_cols,num_nonzeros_diag); hypre_ParCSRBlockMatrixOffd(matrix) = hypre_CSRBlockMatrixCreate(block_size, local_num_rows, num_cols_offd, num_nonzeros_offd); hypre_ParCSRBlockMatrixBlockSize(matrix) = block_size; hypre_ParCSRBlockMatrixGlobalNumRows(matrix) = global_num_rows; hypre_ParCSRBlockMatrixGlobalNumCols(matrix) = global_num_cols; hypre_ParCSRBlockMatrixFirstRowIndex(matrix) = first_row_index; hypre_ParCSRBlockMatrixFirstColDiag(matrix) = first_col_diag; hypre_ParCSRBlockMatrixLastRowIndex(matrix) = first_row_index + local_num_rows - 1; hypre_ParCSRBlockMatrixLastColDiag(matrix) = first_col_diag + local_num_cols - 1; hypre_ParCSRBlockMatrixColMapOffd(matrix) = NULL; hypre_ParCSRBlockMatrixAssumedPartition(matrix) = NULL; /* When NO_GLOBAL_PARTITION is set we could make these null, instead of leaving the range. If that change is made, then when this create is called from functions like the matrix-matrix multiply, be careful not to generate a new partition */ hypre_ParCSRBlockMatrixRowStarts(matrix) = row_starts; hypre_ParCSRBlockMatrixColStarts(matrix) = col_starts; hypre_ParCSRBlockMatrixCommPkg(matrix) = NULL; hypre_ParCSRBlockMatrixCommPkgT(matrix) = NULL; /* set defaults */ hypre_ParCSRBlockMatrixOwnsData(matrix) = 1; hypre_ParCSRBlockMatrixOwnsRowStarts(matrix) = 1; hypre_ParCSRBlockMatrixOwnsColStarts(matrix) = 1; if (row_starts == col_starts) hypre_ParCSRBlockMatrixOwnsColStarts(matrix) = 0; return matrix; }
HYPRE_Int hypre_BoomerAMGBlockCreateNodalA(hypre_ParCSRBlockMatrix *A, HYPRE_Int option, HYPRE_Int diag_option, hypre_ParCSRMatrix **AN_ptr) { MPI_Comm comm = hypre_ParCSRBlockMatrixComm(A); hypre_CSRBlockMatrix *A_diag = hypre_ParCSRBlockMatrixDiag(A); HYPRE_Int *A_diag_i = hypre_CSRBlockMatrixI(A_diag); HYPRE_Real *A_diag_data = hypre_CSRBlockMatrixData(A_diag); HYPRE_Int block_size = hypre_CSRBlockMatrixBlockSize(A_diag); HYPRE_Int bnnz = block_size*block_size; hypre_CSRBlockMatrix *A_offd = hypre_ParCSRMatrixOffd(A); HYPRE_Int *A_offd_i = hypre_CSRBlockMatrixI(A_offd); HYPRE_Real *A_offd_data = hypre_CSRBlockMatrixData(A_offd); HYPRE_Int *A_diag_j = hypre_CSRBlockMatrixJ(A_diag); HYPRE_Int *A_offd_j = hypre_CSRBlockMatrixJ(A_offd); HYPRE_Int *row_starts = hypre_ParCSRBlockMatrixRowStarts(A); HYPRE_Int *col_map_offd = hypre_ParCSRBlockMatrixColMapOffd(A); HYPRE_Int num_nonzeros_diag; HYPRE_Int num_nonzeros_offd = 0; HYPRE_Int num_cols_offd = 0; hypre_ParCSRMatrix *AN; hypre_CSRMatrix *AN_diag; HYPRE_Int *AN_diag_i; HYPRE_Int *AN_diag_j=NULL; HYPRE_Real *AN_diag_data = NULL; hypre_CSRMatrix *AN_offd; HYPRE_Int *AN_offd_i; HYPRE_Int *AN_offd_j = NULL; HYPRE_Real *AN_offd_data = NULL; HYPRE_Int *col_map_offd_AN = NULL; HYPRE_Int *row_starts_AN; hypre_ParCSRCommPkg *comm_pkg = hypre_ParCSRBlockMatrixCommPkg(A); HYPRE_Int num_sends; HYPRE_Int num_recvs; HYPRE_Int *send_procs; HYPRE_Int *send_map_starts; HYPRE_Int *send_map_elmts; HYPRE_Int *recv_procs; HYPRE_Int *recv_vec_starts; hypre_ParCSRCommPkg *comm_pkg_AN = NULL; HYPRE_Int *send_procs_AN = NULL; HYPRE_Int *send_map_starts_AN = NULL; HYPRE_Int *send_map_elmts_AN = NULL; HYPRE_Int *recv_procs_AN = NULL; HYPRE_Int *recv_vec_starts_AN = NULL; HYPRE_Int i; HYPRE_Int ierr = 0; HYPRE_Int num_procs; HYPRE_Int cnt; HYPRE_Int norm_type; HYPRE_Int global_num_nodes; HYPRE_Int num_nodes; HYPRE_Int index, k; HYPRE_Real tmp; HYPRE_Real sum; hypre_MPI_Comm_size(comm,&num_procs); if (!comm_pkg) { hypre_BlockMatvecCommPkgCreate(A); comm_pkg = hypre_ParCSRBlockMatrixCommPkg(A); } norm_type = fabs(option); /* Set up the new matrix AN */ #ifdef HYPRE_NO_GLOBAL_PARTITION row_starts_AN = hypre_CTAlloc(HYPRE_Int, 2); for (i=0; i < 2; i++) { row_starts_AN[i] = row_starts[i]; } #else row_starts_AN = hypre_CTAlloc(HYPRE_Int, num_procs+1); for (i=0; i < num_procs+1; i++) { row_starts_AN[i] = row_starts[i]; } #endif global_num_nodes = hypre_ParCSRBlockMatrixGlobalNumRows(A); num_nodes = hypre_CSRBlockMatrixNumRows(A_diag); /* the diag part */ num_nonzeros_diag = A_diag_i[num_nodes]; AN_diag_i = hypre_CTAlloc(HYPRE_Int, num_nodes+1); for (i=0; i <= num_nodes; i++) { AN_diag_i[i] = A_diag_i[i]; } AN_diag_j = hypre_CTAlloc(HYPRE_Int, num_nonzeros_diag); AN_diag_data = hypre_CTAlloc(HYPRE_Real, num_nonzeros_diag); AN_diag = hypre_CSRMatrixCreate(num_nodes, num_nodes, num_nonzeros_diag); hypre_CSRMatrixI(AN_diag) = AN_diag_i; hypre_CSRMatrixJ(AN_diag) = AN_diag_j; hypre_CSRMatrixData(AN_diag) = AN_diag_data; for (i=0; i< num_nonzeros_diag; i++) { AN_diag_j[i] = A_diag_j[i]; hypre_CSRBlockMatrixBlockNorm(norm_type, &A_diag_data[i*bnnz], &tmp, block_size); AN_diag_data[i] = tmp; } if (diag_option ==1 ) { /* make the diag entry the negative of the sum of off-diag entries (NEED * to get more below!)*/ /* the diagonal is the first element listed in each row - */ for (i=0; i < num_nodes; i++) { index = AN_diag_i[i]; sum = 0.0; for (k = AN_diag_i[i]+1; k < AN_diag_i[i+1]; k++) { sum += AN_diag_data[k]; } AN_diag_data[index] = -sum; } } else if (diag_option == 2) { /* make all diagonal entries negative */ /* the diagonal is the first element listed in each row - */ for (i=0; i < num_nodes; i++) { index = AN_diag_i[i]; AN_diag_data[index] = -AN_diag_data[index]; } } /* copy the commpkg */ if (comm_pkg) { comm_pkg_AN = hypre_CTAlloc(hypre_ParCSRCommPkg,1); hypre_ParCSRCommPkgComm(comm_pkg_AN) = comm; num_sends = hypre_ParCSRCommPkgNumSends(comm_pkg); hypre_ParCSRCommPkgNumSends(comm_pkg_AN) = num_sends; num_recvs = hypre_ParCSRCommPkgNumRecvs(comm_pkg); hypre_ParCSRCommPkgNumRecvs(comm_pkg_AN) = num_recvs; send_procs = hypre_ParCSRCommPkgSendProcs(comm_pkg); send_map_starts = hypre_ParCSRCommPkgSendMapStarts(comm_pkg); send_map_elmts = hypre_ParCSRCommPkgSendMapElmts(comm_pkg); if (num_sends) { send_procs_AN = hypre_CTAlloc(HYPRE_Int, num_sends); send_map_elmts_AN = hypre_CTAlloc(HYPRE_Int, send_map_starts[num_sends]); } send_map_starts_AN = hypre_CTAlloc(HYPRE_Int, num_sends+1); send_map_starts_AN[0] = 0; for (i=0; i < num_sends; i++) { send_procs_AN[i] = send_procs[i]; send_map_starts_AN[i+1] = send_map_starts[i+1]; } cnt = send_map_starts_AN[num_sends]; for (i=0; i< cnt; i++) { send_map_elmts_AN[i] = send_map_elmts[i]; } hypre_ParCSRCommPkgSendProcs(comm_pkg_AN) = send_procs_AN; hypre_ParCSRCommPkgSendMapStarts(comm_pkg_AN) = send_map_starts_AN; hypre_ParCSRCommPkgSendMapElmts(comm_pkg_AN) = send_map_elmts_AN; recv_procs = hypre_ParCSRCommPkgRecvProcs(comm_pkg); recv_vec_starts = hypre_ParCSRCommPkgRecvVecStarts(comm_pkg); recv_vec_starts_AN = hypre_CTAlloc(HYPRE_Int, num_recvs+1); if (num_recvs) recv_procs_AN = hypre_CTAlloc(HYPRE_Int, num_recvs); recv_vec_starts_AN[0] = recv_vec_starts[0]; for (i=0; i < num_recvs; i++) { recv_procs_AN[i] = recv_procs[i]; recv_vec_starts_AN[i+1] = recv_vec_starts[i+1]; } hypre_ParCSRCommPkgRecvProcs(comm_pkg_AN) = recv_procs_AN; hypre_ParCSRCommPkgRecvVecStarts(comm_pkg_AN) = recv_vec_starts_AN; } /* the off-diag part */ num_cols_offd = hypre_CSRBlockMatrixNumCols(A_offd); col_map_offd_AN = hypre_CTAlloc(HYPRE_Int, num_cols_offd); for (i=0; i < num_cols_offd; i++) { col_map_offd_AN[i] = col_map_offd[i]; } num_nonzeros_offd = A_offd_i[num_nodes]; AN_offd_i = hypre_CTAlloc(HYPRE_Int, num_nodes+1); for (i=0; i <= num_nodes; i++) { AN_offd_i[i] = A_offd_i[i]; } AN_offd_j = hypre_CTAlloc(HYPRE_Int, num_nonzeros_offd); AN_offd_data = hypre_CTAlloc(HYPRE_Real, num_nonzeros_offd); for (i=0; i< num_nonzeros_offd; i++) { AN_offd_j[i] = A_offd_j[i]; hypre_CSRBlockMatrixBlockNorm(norm_type, &A_offd_data[i*bnnz], &tmp, block_size); AN_offd_data[i] = tmp; } AN_offd = hypre_CSRMatrixCreate(num_nodes, num_cols_offd, num_nonzeros_offd); hypre_CSRMatrixI(AN_offd) = AN_offd_i; hypre_CSRMatrixJ(AN_offd) = AN_offd_j; hypre_CSRMatrixData(AN_offd) = AN_offd_data; if (diag_option ==1 ) { /* make the diag entry the negative of the sum of off-diag entries (here we are adding the off_diag contribution)*/ /* the diagonal is the first element listed in each row of AN_diag_data - */ for (i=0; i < num_nodes; i++) { sum = 0.0; for (k = AN_offd_i[i]; k < AN_offd_i[i+1]; k++) { sum += AN_offd_data[k]; } index = AN_diag_i[i];/* location of diag entry in data */ AN_diag_data[index] -= sum; /* subtract from current value */ } } /* now create AN */ AN = hypre_ParCSRMatrixCreate(comm, global_num_nodes, global_num_nodes, row_starts_AN, row_starts_AN, num_cols_offd, num_nonzeros_diag, num_nonzeros_offd); /* we already created the diag and offd matrices - so we don't need the ones created above */ hypre_CSRMatrixDestroy(hypre_ParCSRMatrixDiag(AN)); hypre_CSRMatrixDestroy(hypre_ParCSRMatrixOffd(AN)); hypre_ParCSRMatrixDiag(AN) = AN_diag; hypre_ParCSRMatrixOffd(AN) = AN_offd; hypre_ParCSRMatrixColMapOffd(AN) = col_map_offd_AN; hypre_ParCSRMatrixCommPkg(AN) = comm_pkg_AN; *AN_ptr = AN; return (ierr); }