void insert_missing_diags_private(Mat_dh A) { START_FUNC_DH HYPRE_Int *RP = A->rp, *CVAL = A->cval, m = A->m; HYPRE_Int *rp, *cval; double *AVAL = A->aval, *aval; HYPRE_Int i, j, nz = RP[m]+m; HYPRE_Int idx = 0; rp = A->rp = (HYPRE_Int *)MALLOC_DH((1+m)*sizeof(HYPRE_Int)); CHECK_V_ERROR; cval = A->cval = (HYPRE_Int *)MALLOC_DH(nz*sizeof(HYPRE_Int)); CHECK_V_ERROR; aval = A->aval = (double *)MALLOC_DH(nz*sizeof(double)); CHECK_V_ERROR; rp[0] = 0; for (i=0; i<m; ++i) { bool isMissing = true; for (j=RP[i]; j<RP[i+1]; ++j) { cval[idx] = CVAL[j]; aval[idx] = AVAL[j]; ++idx; if (CVAL[j] == i) isMissing = false; } if (isMissing) { cval[idx] = i; aval[idx] = 0.0; ++idx; } rp[i+1] = idx; } FREE_DH(RP); CHECK_V_ERROR; FREE_DH(CVAL); CHECK_V_ERROR; FREE_DH(AVAL); CHECK_V_ERROR; END_FUNC_DH }
void PrintMatUsingGetRow(void* A, HYPRE_Int beg_row, HYPRE_Int m, HYPRE_Int *n2o_row, HYPRE_Int *n2o_col, char *filename) { START_FUNC_DH FILE *fp; HYPRE_Int *o2n_col = NULL, pe, i, j, *cval, len; HYPRE_Int newCol, newRow; double *aval; /* form inverse column permutation */ if (n2o_col != NULL) { o2n_col = (HYPRE_Int*)MALLOC_DH(m*sizeof(HYPRE_Int)); CHECK_V_ERROR; for (i=0; i<m; ++i) o2n_col[n2o_col[i]] = i; } for (pe=0; pe<np_dh; ++pe) { hypre_MPI_Barrier(comm_dh); if (myid_dh == pe) { if (pe == 0) { fp=fopen(filename, "w"); } else { fp=fopen(filename, "a"); } if (fp == NULL) { hypre_sprintf(msgBuf_dh, "can't open %s for writing\n", filename); SET_V_ERROR(msgBuf_dh); } for (i=0; i<m; ++i) { if (n2o_row == NULL) { EuclidGetRow(A, i+beg_row, &len, &cval, &aval); CHECK_V_ERROR; for (j=0; j<len; ++j) { hypre_fprintf(fp, "%i %i %g\n", i+1, cval[j], aval[j]); } EuclidRestoreRow(A, i, &len, &cval, &aval); CHECK_V_ERROR; } else { newRow = n2o_row[i] + beg_row; EuclidGetRow(A, newRow, &len, &cval, &aval); CHECK_V_ERROR; for (j=0; j<len; ++j) { newCol = o2n_col[cval[j]-beg_row] + beg_row; hypre_fprintf(fp, "%i %i %g\n", i+1, newCol, aval[j]); } EuclidRestoreRow(A, i, &len, &cval, &aval); CHECK_V_ERROR; } } fclose(fp); } } if (n2o_col != NULL) { FREE_DH(o2n_col); CHECK_V_ERROR; } END_FUNC_DH }
void convert_triples_to_scr_private(HYPRE_Int m, HYPRE_Int nz, HYPRE_Int *I, HYPRE_Int *J, double *A, HYPRE_Int *rp, HYPRE_Int *cval, double *aval) { START_FUNC_DH HYPRE_Int i; HYPRE_Int *rowCounts; rowCounts = (HYPRE_Int*)MALLOC_DH((m+1)*sizeof(HYPRE_Int)); CHECK_V_ERROR; for (i=0; i<m; ++i) rowCounts[i] = 0; /* count number of entries in each row */ for (i=0; i<nz; ++i) { HYPRE_Int row = I[i]; rowCounts[row] += 1; } /* prefix-sum to form rp[] */ rp[0] = 0; for (i=1; i<=m; ++i) { rp[i] = rp[i-1] + rowCounts[i-1]; } memcpy(rowCounts, rp, (m+1)*sizeof(HYPRE_Int)); /* write SCR arrays */ for (i=0; i<nz; ++i) { HYPRE_Int row = I[i]; HYPRE_Int col = J[i]; double val = A[i]; HYPRE_Int idx = rowCounts[row]; rowCounts[row] += 1; cval[idx] = col; aval[idx] = val; } FREE_DH(rowCounts); CHECK_V_ERROR; END_FUNC_DH }
void iluk_seq(Euclid_dh ctx) { START_FUNC_DH HYPRE_Int *rp, *cval, *diag; HYPRE_Int *CVAL; HYPRE_Int i, j, len, count, col, idx = 0; HYPRE_Int *list, *marker, *fill, *tmpFill; HYPRE_Int temp, m, from = ctx->from, to = ctx->to; HYPRE_Int *n2o_row, *o2n_col, beg_row, beg_rowP; double *AVAL; REAL_DH *work, *aval; Factor_dh F = ctx->F; SubdomainGraph_dh sg = ctx->sg; bool debug = false; if (logFile != NULL && Parser_dhHasSwitch(parser_dh, "-debug_ilu")) debug = true; m = F->m; rp = F->rp; cval = F->cval; fill = F->fill; diag = F->diag; aval = F->aval; work = ctx->work; count = rp[from]; if (sg == NULL) { SET_V_ERROR("subdomain graph is NULL"); } n2o_row = ctx->sg->n2o_row; o2n_col = ctx->sg->o2n_col; beg_row = ctx->sg->beg_row[myid_dh]; beg_rowP = ctx->sg->beg_rowP[myid_dh]; /* allocate and initialize working space */ list = (HYPRE_Int*)MALLOC_DH((m+1)*sizeof(HYPRE_Int)); CHECK_V_ERROR; marker = (HYPRE_Int*)MALLOC_DH(m*sizeof(HYPRE_Int)); CHECK_V_ERROR; tmpFill = (HYPRE_Int*)MALLOC_DH(m*sizeof(HYPRE_Int)); CHECK_V_ERROR; for (i=0; i<m; ++i) marker[i] = -1; /* working space for values */ for (i=0; i<m; ++i) work[i] = 0.0; /* printf_dh("====================== starting iluk_seq; level= %i\n\n", ctx->level); */ /*---------- main loop ----------*/ for (i=from; i<to; ++i) { HYPRE_Int row = n2o_row[i]; /* local row number */ HYPRE_Int globalRow = row+beg_row; /* global row number */ /*hypre_fprintf(logFile, "--------------------------------- localRow= %i\n", 1+i); */ if (debug) { hypre_fprintf(logFile, "ILU_seq ================================= starting local row: %i, (global= %i) level= %i\n", i+1, i+1+sg->beg_rowP[myid_dh], ctx->level); } EuclidGetRow(ctx->A, globalRow, &len, &CVAL, &AVAL); CHECK_V_ERROR; /* compute scaling value for row(i) */ if (ctx->isScaled) { compute_scaling_private(i, len, AVAL, ctx); CHECK_V_ERROR; } /* Compute symbolic factor for row(i); this also performs sparsification */ count = symbolic_row_private(i, list, marker, tmpFill, len, CVAL, AVAL, o2n_col, ctx, debug); CHECK_V_ERROR; /* Ensure adequate storage; reallocate, if necessary. */ if (idx + count > F->alloc) { Factor_dhReallocate(F, idx, count); CHECK_V_ERROR; SET_INFO("REALLOCATED from ilu_seq"); cval = F->cval; fill = F->fill; aval = F->aval; } /* Copy factored symbolic row to permanent storage */ col = list[m]; while (count--) { cval[idx] = col; fill[idx] = tmpFill[col]; ++idx; /*hypre_fprintf(logFile, " col= %i\n", 1+col); */ col = list[col]; } /* add row-pointer to start of next row. */ rp[i+1] = idx; /* Insert pointer to diagonal */ temp = rp[i]; while (cval[temp] != i) ++temp; diag[i] = temp; /*hypre_fprintf(logFile, " diag[i]= %i\n", diag); */ /* compute numeric factor for current row */ numeric_row_private(i, len, CVAL, AVAL, work, o2n_col, ctx, debug); CHECK_V_ERROR EuclidRestoreRow(ctx->A, globalRow, &len, &CVAL, &AVAL); CHECK_V_ERROR; /* Copy factored numeric row to permanent storage, and re-zero work vector */ if (debug) { hypre_fprintf(logFile, "ILU_seq: "); for (j=rp[i]; j<rp[i+1]; ++j) { col = cval[j]; aval[j] = work[col]; work[col] = 0.0; hypre_fprintf(logFile, "%i,%i,%g ; ", 1+cval[j], fill[j], aval[j]); fflush(logFile); } hypre_fprintf(logFile, "\n"); } else { for (j=rp[i]; j<rp[i+1]; ++j) { col = cval[j]; aval[j] = work[col]; work[col] = 0.0; } } /* check for zero diagonal */ if (! aval[diag[i]]) { hypre_sprintf(msgBuf_dh, "zero diagonal in local row %i", i+1); SET_V_ERROR(msgBuf_dh); } } FREE_DH(list); CHECK_V_ERROR; FREE_DH(tmpFill); CHECK_V_ERROR; FREE_DH(marker); CHECK_V_ERROR; /* adjust column indices back to global */ if (beg_rowP) { HYPRE_Int start = rp[from]; HYPRE_Int stop = rp[to]; for (i=start; i<stop; ++i) cval[i] += beg_rowP; } /* for debugging: this is so the Print methods will work, even if F hasn't been fully factored */ for (i=to+1; i<m; ++i) rp[i] = 0; END_FUNC_DH }
void ilut_seq(Euclid_dh ctx) { START_FUNC_DH HYPRE_Int *rp, *cval, *diag, *CVAL; HYPRE_Int i, len, count, col, idx = 0; HYPRE_Int *list, *marker; HYPRE_Int temp, m, from, to; HYPRE_Int *n2o_row, *o2n_col, beg_row, beg_rowP; double *AVAL, droptol; REAL_DH *work, *aval, val; Factor_dh F = ctx->F; SubdomainGraph_dh sg = ctx->sg; bool debug = false; if (logFile != NULL && Parser_dhHasSwitch(parser_dh, "-debug_ilu")) debug = true; m = F->m; rp = F->rp; cval = F->cval; diag = F->diag; aval = F->aval; work = ctx->work; from = ctx->from; to = ctx->to; count = rp[from]; droptol = ctx->droptol; if (sg == NULL) { SET_V_ERROR("subdomain graph is NULL"); } n2o_row = ctx->sg->n2o_row; o2n_col = ctx->sg->o2n_col; beg_row = ctx->sg->beg_row[myid_dh]; beg_rowP = ctx->sg->beg_rowP[myid_dh]; /* allocate and initialize working space */ list = (HYPRE_Int*)MALLOC_DH((m+1)*sizeof(HYPRE_Int)); CHECK_V_ERROR; marker = (HYPRE_Int*)MALLOC_DH(m*sizeof(HYPRE_Int)); CHECK_V_ERROR; for (i=0; i<m; ++i) marker[i] = -1; rp[0] = 0; /* working space for values */ for (i=0; i<m; ++i) work[i] = 0.0; /* ----- main loop start ----- */ for (i=from; i<to; ++i) { HYPRE_Int row = n2o_row[i]; /* local row number */ HYPRE_Int globalRow = row + beg_row; /* global row number */ EuclidGetRow(ctx->A, globalRow, &len, &CVAL, &AVAL); CHECK_V_ERROR; /* compute scaling value for row(i) */ compute_scaling_private(i, len, AVAL, ctx); CHECK_V_ERROR; /* compute factor for row i */ count = ilut_row_private(i, list, o2n_col, marker, len, CVAL, AVAL, work, ctx, debug); CHECK_V_ERROR; EuclidRestoreRow(ctx->A, globalRow, &len, &CVAL, &AVAL); CHECK_V_ERROR; /* Ensure adequate storage; reallocate, if necessary. */ if (idx + count > F->alloc) { Factor_dhReallocate(F, idx, count); CHECK_V_ERROR; SET_INFO("REALLOCATED from ilu_seq"); cval = F->cval; aval = F->aval; } /* Copy factored row to permanent storage, apply 2nd drop test, and re-zero work vector */ col = list[m]; while (count--) { val = work[col]; if (col == i || fabs(val) > droptol) { cval[idx] = col; aval[idx++] = val; work[col] = 0.0; } col = list[col]; } /* add row-pointer to start of next row. */ rp[i+1] = idx; /* Insert pointer to diagonal */ temp = rp[i]; while (cval[temp] != i) ++temp; diag[i] = temp; /* check for zero diagonal */ if (! aval[diag[i]]) { hypre_sprintf(msgBuf_dh, "zero diagonal in local row %i", i+1); SET_V_ERROR(msgBuf_dh); } } /* --------- main loop end --------- */ /* adjust column indices back to global */ if (beg_rowP) { HYPRE_Int start = rp[from]; HYPRE_Int stop = rp[to]; for (i=start; i<stop; ++i) cval[i] += beg_rowP; } FREE_DH(list); FREE_DH(marker); END_FUNC_DH }
void iluk_seq_block(Euclid_dh ctx) { START_FUNC_DH HYPRE_Int *rp, *cval, *diag; HYPRE_Int *CVAL; HYPRE_Int h, i, j, len, count, col, idx = 0; HYPRE_Int *list, *marker, *fill, *tmpFill; HYPRE_Int temp, m; HYPRE_Int *n2o_row, *o2n_col, *beg_rowP, *n2o_sub, blocks; HYPRE_Int *row_count, *dummy = NULL, dummy2[1]; double *AVAL; REAL_DH *work, *aval; Factor_dh F = ctx->F; SubdomainGraph_dh sg = ctx->sg; bool bj = false, constrained = false; HYPRE_Int discard = 0; HYPRE_Int gr = -1; /* globalRow */ bool debug = false; if (logFile != NULL && Parser_dhHasSwitch(parser_dh, "-debug_ilu")) debug = true; /*hypre_fprintf(stderr, "====================== starting iluk_seq_block; level= %i\n\n", ctx->level); */ if (!strcmp(ctx->algo_par, "bj")) bj = true; constrained = ! Parser_dhHasSwitch(parser_dh, "-unconstrained"); m = F->m; rp = F->rp; cval = F->cval; fill = F->fill; diag = F->diag; aval = F->aval; work = ctx->work; if (sg != NULL) { n2o_row = sg->n2o_row; o2n_col = sg->o2n_col; row_count = sg->row_count; /* beg_row = sg->beg_row ; */ beg_rowP = sg->beg_rowP; n2o_sub = sg->n2o_sub; blocks = sg->blocks; } else { dummy = (HYPRE_Int*)MALLOC_DH(m*sizeof(HYPRE_Int)); CHECK_V_ERROR; for (i=0; i<m; ++i) dummy[i] = i; n2o_row = dummy; o2n_col = dummy; dummy2[0] = m; row_count = dummy2; /* beg_row = 0; */ beg_rowP = dummy; n2o_sub = dummy; blocks = 1; } /* allocate and initialize working space */ list = (HYPRE_Int*)MALLOC_DH((m+1)*sizeof(HYPRE_Int)); CHECK_V_ERROR; marker = (HYPRE_Int*)MALLOC_DH(m*sizeof(HYPRE_Int)); CHECK_V_ERROR; tmpFill = (HYPRE_Int*)MALLOC_DH(m*sizeof(HYPRE_Int)); CHECK_V_ERROR; for (i=0; i<m; ++i) marker[i] = -1; /* working space for values */ for (i=0; i<m; ++i) work[i] = 0.0; /*---------- main loop ----------*/ for (h=0; h<blocks; ++h) { /* 1st and last row in current block, with respect to A */ HYPRE_Int curBlock = n2o_sub[h]; HYPRE_Int first_row = beg_rowP[curBlock]; HYPRE_Int end_row = first_row + row_count[curBlock]; if (debug) { hypre_fprintf(logFile, "\n\nILU_seq BLOCK: %i @@@@@@@@@@@@@@@ \n", curBlock); } for (i=first_row; i<end_row; ++i) { HYPRE_Int row = n2o_row[i]; ++gr; if (debug) { hypre_fprintf(logFile, "ILU_seq global: %i local: %i =================================\n", 1+gr, 1+i-first_row); } /*prinft("first_row= %i end_row= %i\n", first_row, end_row); */ EuclidGetRow(ctx->A, row, &len, &CVAL, &AVAL); CHECK_V_ERROR; /* compute scaling value for row(i) */ if (ctx->isScaled) { compute_scaling_private(i, len, AVAL, ctx); CHECK_V_ERROR; } /* Compute symbolic factor for row(i); this also performs sparsification */ count = symbolic_row_private(i, list, marker, tmpFill, len, CVAL, AVAL, o2n_col, ctx, debug); CHECK_V_ERROR; /* Ensure adequate storage; reallocate, if necessary. */ if (idx + count > F->alloc) { Factor_dhReallocate(F, idx, count); CHECK_V_ERROR; SET_INFO("REALLOCATED from ilu_seq"); cval = F->cval; fill = F->fill; aval = F->aval; } /* Copy factored symbolic row to permanent storage */ col = list[m]; while (count--) { /* constrained pilu */ if (constrained && !bj) { if (col >= first_row && col < end_row) { cval[idx] = col; fill[idx] = tmpFill[col]; ++idx; } else { if (check_constraint_private(ctx, curBlock, col)) { cval[idx] = col; fill[idx] = tmpFill[col]; ++idx; } else { ++discard; } } col = list[col]; } /* block jacobi case */ else if (bj) { if (col >= first_row && col < end_row) { cval[idx] = col; fill[idx] = tmpFill[col]; ++idx; } else { ++discard; } col = list[col]; } /* general case */ else { cval[idx] = col; fill[idx] = tmpFill[col]; ++idx; col = list[col]; } } /* add row-pointer to start of next row. */ rp[i+1] = idx; /* Insert pointer to diagonal */ temp = rp[i]; while (cval[temp] != i) ++temp; diag[i] = temp; /* compute numeric factor for current row */ numeric_row_private(i, len, CVAL, AVAL, work, o2n_col, ctx, debug); CHECK_V_ERROR EuclidRestoreRow(ctx->A, row, &len, &CVAL, &AVAL); CHECK_V_ERROR; /* Copy factored numeric row to permanent storage, and re-zero work vector */ if (debug) { hypre_fprintf(logFile, "ILU_seq: "); for (j=rp[i]; j<rp[i+1]; ++j) { col = cval[j]; aval[j] = work[col]; work[col] = 0.0; hypre_fprintf(logFile, "%i,%i,%g ; ", 1+cval[j], fill[j], aval[j]); } hypre_fprintf(logFile, "\n"); } /* normal operation */ else { for (j=rp[i]; j<rp[i+1]; ++j) { col = cval[j]; aval[j] = work[col]; work[col] = 0.0; } } /* check for zero diagonal */ if (! aval[diag[i]]) { hypre_sprintf(msgBuf_dh, "zero diagonal in local row %i", i+1); SET_V_ERROR(msgBuf_dh); } } } /* hypre_printf("bj= %i constrained= %i discarded= %i\n", bj, constrained, discard); */ if (dummy != NULL) { FREE_DH(dummy); CHECK_V_ERROR; } FREE_DH(list); CHECK_V_ERROR; FREE_DH(tmpFill); CHECK_V_ERROR; FREE_DH(marker); CHECK_V_ERROR; END_FUNC_DH }
void iluk_mpi_bj(Euclid_dh ctx) { START_FUNC_DH HYPRE_Int *rp, *cval, *diag; HYPRE_Int *CVAL; HYPRE_Int i, j, len, count, col, idx = 0; HYPRE_Int *list, *marker, *fill, *tmpFill; HYPRE_Int temp, m, from = ctx->from, to = ctx->to; HYPRE_Int *n2o_row, *o2n_col; HYPRE_Int first_row, last_row; double *AVAL; REAL_DH *work, *aval; Factor_dh F = ctx->F; SubdomainGraph_dh sg = ctx->sg; if (ctx->F == NULL) { SET_V_ERROR("ctx->F is NULL"); } if (ctx->F->rp == NULL) { SET_V_ERROR("ctx->F->rp is NULL"); } /* printf_dh("====================== starting iluk_mpi_bj; level= %i\n\n", ctx->level); */ m = F->m; rp = F->rp; cval = F->cval; fill = F->fill; diag = F->diag; aval = F->aval; work = ctx->work; n2o_row = sg->n2o_row; o2n_col = sg->o2n_col; /* allocate and initialize working space */ list = (HYPRE_Int*)MALLOC_DH((m+1)*sizeof(HYPRE_Int)); CHECK_V_ERROR; marker = (HYPRE_Int*)MALLOC_DH(m*sizeof(HYPRE_Int)); CHECK_V_ERROR; tmpFill = (HYPRE_Int*)MALLOC_DH(m*sizeof(HYPRE_Int)); CHECK_V_ERROR; for (i=0; i<m; ++i) { marker[i] = -1; work[i] = 0.0; } /*---------- main loop ----------*/ /* global numbers of first and last locally owned rows, with respect to A */ first_row = sg->beg_row[myid_dh]; last_row = first_row + sg->row_count[myid_dh]; for (i=from; i<to; ++i) { HYPRE_Int row = n2o_row[i]; /* local row number */ HYPRE_Int globalRow = row + first_row; /* global row number */ EuclidGetRow(ctx->A, globalRow, &len, &CVAL, &AVAL); CHECK_V_ERROR; /* compute scaling value for row(i) */ if (ctx->isScaled) { compute_scaling_private(i, len, AVAL, ctx); CHECK_V_ERROR; } /* Compute symbolic factor for row(i); this also performs sparsification */ count = symbolic_row_private(i, first_row, last_row, list, marker, tmpFill, len, CVAL, AVAL, o2n_col, ctx); CHECK_V_ERROR; /* Ensure adequate storage; reallocate, if necessary. */ if (idx + count > F->alloc) { Factor_dhReallocate(F, idx, count); CHECK_V_ERROR; SET_INFO("REALLOCATED from lu_mpi_bj"); cval = F->cval; fill = F->fill; aval = F->aval; } /* Copy factored symbolic row to permanent storage */ col = list[m]; while (count--) { cval[idx] = col; fill[idx] = tmpFill[col]; ++idx; col = list[col]; } /* add row-pointer to start of next row. */ rp[i+1] = idx; /* Insert pointer to diagonal */ temp = rp[i]; while (cval[temp] != i) ++temp; diag[i] = temp; /* compute numeric factor for current row */ numeric_row_private(i, first_row, last_row, len, CVAL, AVAL, work, o2n_col, ctx); CHECK_V_ERROR EuclidRestoreRow(ctx->A, globalRow, &len, &CVAL, &AVAL); CHECK_V_ERROR; /* Copy factored numeric row to permanent storage, and re-zero work vector */ for (j=rp[i]; j<rp[i+1]; ++j) { col = cval[j]; aval[j] = work[col]; work[col] = 0.0; } /* check for zero diagonal */ if (! aval[diag[i]]) { hypre_sprintf(msgBuf_dh, "zero diagonal in local row %i", i+1); SET_V_ERROR(msgBuf_dh); } } FREE_DH(list); CHECK_V_ERROR; FREE_DH(tmpFill); CHECK_V_ERROR; FREE_DH(marker); CHECK_V_ERROR; END_FUNC_DH }
void bicgstab_euclid(Mat_dh A, Euclid_dh ctx, double *x, double *b, HYPRE_Int *itsOUT) { START_FUNC_DH HYPRE_Int its, m = ctx->m; bool monitor; HYPRE_Int maxIts = ctx->maxIts; double atol = ctx->atol, rtol = ctx->rtol; /* scalars */ double alpha, alpha_1, beta_1, widget, widget_1, rho_1, rho_2, s_norm, eps, exit_a, b_iprod, r_iprod; /* vectors */ double *t, *s, *s_hat, *v, *p, *p_hat, *r, *r_hat; monitor = Parser_dhHasSwitch(parser_dh, "-monitor"); /* allocate working space */ t = (double*)MALLOC_DH(m*sizeof(double)); s = (double*)MALLOC_DH(m*sizeof(double)); s_hat = (double*)MALLOC_DH(m*sizeof(double)); v = (double*)MALLOC_DH(m*sizeof(double)); p = (double*)MALLOC_DH(m*sizeof(double)); p_hat = (double*)MALLOC_DH(m*sizeof(double)); r = (double*)MALLOC_DH(m*sizeof(double)); r_hat = (double*)MALLOC_DH(m*sizeof(double)); /* r = b - Ax */ Mat_dhMatVec(A, x, s); /* s = Ax */ CopyVec(m, b, r); /* r = b */ Axpy(m, -1.0, s, r); /* r = b-Ax */ CopyVec(m, r, r_hat); /* r_hat = r */ /* compute stopping criteria */ b_iprod = InnerProd(m, b, b); CHECK_V_ERROR; exit_a = atol*atol*b_iprod; CHECK_V_ERROR; /* absolute stopping criteria */ eps = rtol*rtol*b_iprod; /* relative stoping criteria (residual reduction) */ its = 0; while(1) { ++its; rho_1 = InnerProd(m, r_hat, r); if (rho_1 == 0) { SET_V_ERROR("(r_hat . r) = 0; method fails"); } if (its == 1) { CopyVec(m, r, p); /* p = r_0 */ CHECK_V_ERROR; } else { beta_1 = (rho_1/rho_2)*(alpha_1/widget_1); /* p_i = r_(i-1) + beta_(i-1)*( p_(i-1) - w_(i-1)*v_(i-1) ) */ Axpy(m, -widget_1, v, p); CHECK_V_ERROR; ScaleVec(m, beta_1, p); CHECK_V_ERROR; Axpy(m, 1.0, r, p); CHECK_V_ERROR; } /* solve M*p_hat = p_i */ Euclid_dhApply(ctx, p, p_hat); CHECK_V_ERROR; /* v_i = A*p_hat */ Mat_dhMatVec(A, p_hat, v); CHECK_V_ERROR; /* alpha_i = rho_(i-1) / (r_hat^T . v_i ) */ { double tmp = InnerProd(m, r_hat, v); CHECK_V_ERROR; alpha = rho_1/tmp; } /* s = r_(i-1) - alpha_i*v_i */ CopyVec(m, r, s); CHECK_V_ERROR; Axpy(m, -alpha, v, s); CHECK_V_ERROR; /* check norm of s; if small enough: * set x_i = x_(i-1) + alpha_i*p_i and stop. * (Actually, we use the square of the norm) */ s_norm = InnerProd(m, s, s); if (s_norm < exit_a) { SET_INFO("reached absolute stopping criteria"); break; } /* solve M*s_hat = s */ Euclid_dhApply(ctx, s, s_hat); CHECK_V_ERROR; /* t = A*s_hat */ Mat_dhMatVec(A, s_hat, t); CHECK_V_ERROR; /* w_i = (t . s)/(t . t) */ { double tmp1, tmp2; tmp1 = InnerProd(m, t, s); CHECK_V_ERROR; tmp2 = InnerProd(m, t, t); CHECK_V_ERROR; widget = tmp1/tmp2; } /* x_i = x_(i-1) + alpha_i*p_hat + w_i*s_hat */ Axpy(m, alpha, p_hat, x); CHECK_V_ERROR; Axpy(m, widget, s_hat, x); CHECK_V_ERROR; /* r_i = s - w_i*t */ CopyVec(m, s, r); CHECK_V_ERROR; Axpy(m, -widget, t, r); CHECK_V_ERROR; /* check convergence; continue if necessary; * for continuation it is necessary thea w != 0. */ r_iprod = InnerProd(m, r, r); CHECK_V_ERROR; if (r_iprod < eps) { SET_INFO("stipulated residual reduction achieved"); break; } /* monitor convergence */ if (monitor && myid_dh == 0) { hypre_fprintf(stderr, "[it = %i] %e\n", its, sqrt(r_iprod/b_iprod)); } /* prepare for next iteration */ rho_2 = rho_1; widget_1 = widget; alpha_1 = alpha; if (its >= maxIts) { its = -its; break; } } *itsOUT = its; FREE_DH(t); FREE_DH(s); FREE_DH(s_hat); FREE_DH(v); FREE_DH(p); FREE_DH(p_hat); FREE_DH(r); FREE_DH(r_hat); END_FUNC_DH }
void cg_euclid(Mat_dh A, Euclid_dh ctx, double *x, double *b, HYPRE_Int *itsOUT) { START_FUNC_DH HYPRE_Int its, m = A->m; double *p, *r, *s; double alpha, beta, gamma, gamma_old, eps, bi_prod, i_prod; bool monitor; HYPRE_Int maxIts = ctx->maxIts; /* double atol = ctx->atol */ double rtol = ctx->rtol; monitor = Parser_dhHasSwitch(parser_dh, "-monitor"); /* compute square of absolute stopping threshold */ /* bi_prod = <b,b> */ bi_prod = InnerProd(m, b, b); CHECK_V_ERROR; eps = (rtol*rtol)*bi_prod; p = (double *) MALLOC_DH(m * sizeof(double)); s = (double *) MALLOC_DH(m * sizeof(double)); r = (double *) MALLOC_DH(m * sizeof(double)); /* r = b - Ax */ Mat_dhMatVec(A, x, r); /* r = Ax */ CHECK_V_ERROR; ScaleVec(m, -1.0, r); /* r = b */ CHECK_V_ERROR; Axpy(m, 1.0, b, r); /* r = r + b */ CHECK_V_ERROR; /* solve Mp = r */ Euclid_dhApply(ctx, r, p); CHECK_V_ERROR; /* gamma = <r,p> */ gamma = InnerProd(m, r, p); CHECK_V_ERROR; its = 0; while (1) { ++its; /* s = A*p */ Mat_dhMatVec(A, p, s); CHECK_V_ERROR; /* alpha = gamma / <s,p> */ { double tmp = InnerProd(m, s, p); CHECK_V_ERROR; alpha = gamma / tmp; gamma_old = gamma; } /* x = x + alpha*p */ Axpy(m, alpha, p, x); CHECK_V_ERROR; /* r = r - alpha*s */ Axpy(m, -alpha, s, r); CHECK_V_ERROR; /* solve Ms = r */ Euclid_dhApply(ctx, r, s); CHECK_V_ERROR; /* gamma = <r,s> */ gamma = InnerProd(m, r, s); CHECK_V_ERROR; /* set i_prod for convergence test */ i_prod = InnerProd(m, r, r); CHECK_V_ERROR; if (monitor && myid_dh == 0) { hypre_fprintf(stderr, "iter = %i rel. resid. norm: %e\n", its, sqrt(i_prod/bi_prod)); } /* check for convergence */ if (i_prod < eps) break; /* beta = gamma / gamma_old */ beta = gamma / gamma_old; /* p = s + beta p */ ScaleVec(m, beta, p); CHECK_V_ERROR; Axpy(m, 1.0, s, p); CHECK_V_ERROR; if (its >= maxIts) { its = -its; break; } } *itsOUT = its; FREE_DH(p); FREE_DH(s); FREE_DH(r); END_FUNC_DH }
void ExternalRows_dhDestroy (ExternalRows_dh er) { START_FUNC_DH int i; for (i = 0; i < MAX_MPI_TASKS; ++i) { if (er->rcv_row_lengths[i] != NULL) { FREE_DH (er->rcv_row_lengths[i]); CHECK_V_ERROR; } if (er->rcv_row_numbers[i] != NULL) { FREE_DH (er->rcv_row_numbers[i]); CHECK_V_ERROR; } } if (er->cvalExt != NULL) { FREE_DH (er->cvalExt); CHECK_V_ERROR; } if (er->fillExt != NULL) { FREE_DH (er->fillExt); CHECK_V_ERROR; } if (er->avalExt != NULL) { FREE_DH (er->avalExt); CHECK_V_ERROR; } if (er->my_row_counts != NULL) { FREE_DH (er->my_row_counts); CHECK_V_ERROR; } if (er->my_row_numbers != NULL) { FREE_DH (er->my_row_numbers); CHECK_V_ERROR; } if (er->cvalSend != NULL) { FREE_DH (er->cvalSend); CHECK_V_ERROR; } if (er->fillSend != NULL) { FREE_DH (er->fillSend); CHECK_V_ERROR; } if (er->avalSend != NULL) { FREE_DH (er->avalSend); CHECK_V_ERROR; } if (er->rowLookup != NULL) { Hash_dhDestroy (er->rowLookup); CHECK_V_ERROR; } FREE_DH (er); CHECK_V_ERROR; END_FUNC_DH}
void MatGenFD_Destroy(MatGenFD mg) { START_FUNC_DH FREE_DH(mg); CHECK_V_ERROR; END_FUNC_DH }
void mat_dh_read_triples_private(HYPRE_Int ignore, HYPRE_Int *mOUT, HYPRE_Int **rpOUT, HYPRE_Int **cvalOUT, double **avalOUT, FILE* fp) { START_FUNC_DH HYPRE_Int m, n, nz, items, i, j; HYPRE_Int idx = 0; HYPRE_Int *cval, *rp, *I, *J; double *aval, *A, v; char junk[MAX_JUNK]; fpos_t fpos; /* skip over header */ if (ignore && myid_dh == 0) { hypre_printf("mat_dh_read_triples_private:: ignoring following header lines:\n"); hypre_printf("--------------------------------------------------------------\n"); for (i=0; i<ignore; ++i) { fgets(junk, MAX_JUNK, fp); hypre_printf("%s", junk); } hypre_printf("--------------------------------------------------------------\n"); if (fgetpos(fp, &fpos)) SET_V_ERROR("fgetpos failed!"); hypre_printf("\nmat_dh_read_triples_private::1st two non-ignored lines:\n"); hypre_printf("--------------------------------------------------------------\n"); for (i=0; i<2; ++i) { fgets(junk, MAX_JUNK, fp); hypre_printf("%s", junk); } hypre_printf("--------------------------------------------------------------\n"); if (fsetpos(fp, &fpos)) SET_V_ERROR("fsetpos failed!"); } if (feof(fp)) hypre_printf("trouble!"); /* determine matrix dimensions */ m=n=nz=0; while (!feof(fp)) { items = hypre_fscanf(fp,"%d %d %lg",&i,&j,&v); if (items != 3) { break; } ++nz; if (i > m) m = i; if (j > n) n = j; } if (myid_dh == 0) { hypre_printf("mat_dh_read_triples_private: m= %i nz= %i\n", m, nz); } /* reset file, and skip over header again */ rewind(fp); for (i=0; i<ignore; ++i) { fgets(junk, MAX_JUNK, fp); } /* error check for squareness */ if (m != n) { hypre_sprintf(msgBuf_dh, "matrix is not square; row= %i, cols= %i", m, n); SET_V_ERROR(msgBuf_dh); } *mOUT = m; /* allocate storage */ rp = *rpOUT = (HYPRE_Int*)MALLOC_DH((m+1)*sizeof(HYPRE_Int)); CHECK_V_ERROR; cval = *cvalOUT = (HYPRE_Int*)MALLOC_DH(nz*sizeof(HYPRE_Int)); CHECK_V_ERROR; aval = *avalOUT = (double*)MALLOC_DH(nz*sizeof(double)); CHECK_V_ERROR; I = (HYPRE_Int*)MALLOC_DH(nz*sizeof(HYPRE_Int)); CHECK_V_ERROR; J = (HYPRE_Int*)MALLOC_DH(nz*sizeof(HYPRE_Int)); CHECK_V_ERROR; A = (double*)MALLOC_DH(nz*sizeof(double)); CHECK_V_ERROR; /* read <row, col, value> triples into arrays */ while (!feof(fp)) { items = hypre_fscanf(fp,"%d %d %lg",&i,&j,&v); if (items < 3) break; j--; i--; I[idx] = i; J[idx] = j; A[idx] = v; ++idx; } /* convert from triples to sparse-compressed-row storage */ convert_triples_to_scr_private(m, nz, I, J, A, rp, cval, aval); CHECK_V_ERROR; /* if matrix is triangular */ { HYPRE_Int type; type = isTriangular(m, rp, cval); CHECK_V_ERROR; if (type == IS_UPPER_TRI) { hypre_printf("CAUTION: matrix is upper triangular; converting to full\n"); } else if (type == IS_LOWER_TRI) { hypre_printf("CAUTION: matrix is lower triangular; converting to full\n"); } if (type == IS_UPPER_TRI || type == IS_LOWER_TRI) { make_full_private(m, &rp, &cval, &aval); CHECK_V_ERROR; } } *rpOUT = rp; *cvalOUT = cval; *avalOUT = aval; FREE_DH(I); CHECK_V_ERROR; FREE_DH(J); CHECK_V_ERROR; FREE_DH(A); CHECK_V_ERROR; END_FUNC_DH }
void destroy_nat_ordering_private(HYPRE_Int *p) { START_FUNC_DH FREE_DH(p); CHECK_V_ERROR; END_FUNC_DH }
void mat_dh_print_graph_private(HYPRE_Int m, HYPRE_Int beg_row, HYPRE_Int *rp, HYPRE_Int *cval, double *aval, HYPRE_Int *n2o, HYPRE_Int *o2n, Hash_i_dh hash, FILE* fp) { START_FUNC_DH HYPRE_Int i, j, row, col; bool private_n2o = false; bool private_hash = false; HYPRE_Int *work = NULL; work = (HYPRE_Int*)MALLOC_DH(m*sizeof(HYPRE_Int)); CHECK_V_ERROR; if (n2o == NULL) { private_n2o = true; create_nat_ordering_private(m, &n2o); CHECK_V_ERROR; create_nat_ordering_private(m, &o2n); CHECK_V_ERROR; } if (hash == NULL) { private_hash = true; Hash_i_dhCreate(&hash, -1); CHECK_V_ERROR; } for (i=0; i<m; ++i) { for (j=0; j<m; ++j) work[j] = 0; row = n2o[i]; for (j=rp[row]; j<rp[row+1]; ++j) { col = cval[j]; /* local column */ if (col >= beg_row || col < beg_row+m) { col = o2n[col]; } /* nonlocal column: get permutation from hash table */ else { HYPRE_Int tmp = col; tmp = Hash_i_dhLookup(hash, col); CHECK_V_ERROR; if (tmp == -1) { hypre_sprintf(msgBuf_dh, "beg_row= %i m= %i; nonlocal column= %i not in hash table", beg_row, m, col); SET_V_ERROR(msgBuf_dh); } else { col = tmp; } } work[col] = 1; } for (j=0; j<m; ++j) { if (work[j]) { hypre_fprintf(fp, " x "); } else { hypre_fprintf(fp, " "); } } hypre_fprintf(fp, "\n"); } if (private_n2o) { destroy_nat_ordering_private(n2o); CHECK_V_ERROR; destroy_nat_ordering_private(o2n); CHECK_V_ERROR; } if (private_hash) { Hash_i_dhDestroy(hash); CHECK_V_ERROR; } if (work != NULL) { FREE_DH(work); CHECK_V_ERROR; } END_FUNC_DH }