Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

ILU as a preconditioner and as a stand-alone solver #52

Draft
wants to merge 2 commits into
base: master
Choose a base branch
from
Draft
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
89 changes: 89 additions & 0 deletions src/HypreSystem.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -57,6 +57,9 @@ void HypreSystem::setup_precon_and_solver() {

if (preconditioner == "boomeramg") {
setup_boomeramg_precond();
}
else if (preconditioner == "ilu"){
setup_ilu_precond();
} else if (preconditioner == "none") {
usePrecond_ = false;
} else {
Expand All @@ -75,6 +78,8 @@ void HypreSystem::setup_precon_and_solver() {
setup_boomeramg_solver();
} else if (!method.compare("cogmres")) {
setup_cogmres();
} else if (!method.compare("ilu")) {
setup_ilu();
} else {
throw std::runtime_error("Invalid option for solver method provided: " +
method);
Expand Down Expand Up @@ -286,6 +291,49 @@ void HypreSystem::setup_boomeramg_precond() {
precondDestroyPtr_ = &HYPRE_BoomerAMGDestroy;
}

void HypreSystem::setup_ilu_precond() {
YAML::Node node = inpfile_["ilu_preconditioner_settings"];

HYPRE_ILUCreate(&precond_);
HYPRE_ILUSetType(precond_, get_optional(node, "ilu_type", 0));
HYPRE_ILUSetMaxIter(precond_, get_optional(node, "max_iterations", 1));
HYPRE_ILUSetTol(precond_, get_optional(node, "tolerance", 0.0));
HYPRE_ILUSetLocalReordering(precond_, get_optional(node, "local_reordering", 0));
HYPRE_ILUSetPrintLevel(precond_, get_optional(node, "print_level", 1));

// ILUK parameters
HYPRE_ILUSetLevelOfFill(precond_, get_optional(node, "fill", 0));

// ILUT parameters
HYPRE_ILUSetMaxNnzPerRow(precond_,get_optional(node, "max_nnz_per_row", 1000));
HYPRE_ILUSetDropThreshold(precond_,get_optional(node, "drop_threshold", 1.0e-2));

// 0 : Non-iterative algorithm (default)
// 1 : Asynchronous with in-place storage
// 2 : Asynchronous with explicit storage splitting
// 3 : Synchronous with explicit storage splitting
// 4 : Semi-synchronous with explicit storage splitting
// Iterative ILU is available only for zero fill-in and it depends on rocSparse
HYPRE_ILUSetIterativeSetupType(precond_,
get_optional(node, "algorithm_type", 0));
HYPRE_ILUSetIterativeSetupMaxIter(
precond_, get_optional(node, "max_ilu_iterations", 1));
HYPRE_ILUSetIterativeSetupTolerance(
precond_, get_optional(node, "iterative_ilu_tolerance", 1e-5));

// 0: iterative
// 1: direct (default)
HYPRE_ILUSetTriSolve(precond_, get_optional(node, "trisolve", 1));

// Jacobi iterations for lower and upper triangular solves
HYPRE_ILUSetLowerJacobiIters(precond_,get_optional(node, "lower_jacobi_iters", 5));
HYPRE_ILUSetUpperJacobiIters(precond_,get_optional(node, "upper_jacobi_iters", 5));

precondSetupPtr_ = &HYPRE_ILUSetup;
precondSolvePtr_ = &HYPRE_ILUSolve;
precondDestroyPtr_ = &HYPRE_ILUDestroy;
}

void HypreSystem::setup_cogmres() {
YAML::Node node = inpfile_["solver_settings"];

Expand Down Expand Up @@ -371,6 +419,47 @@ void HypreSystem::setup_cg() {
solverSolvePtr_ = &HYPRE_ParCSRPCGSolve;
}

void HypreSystem::setup_ilu() {
YAML::Node node = inpfile_["solver_settings"];
HYPRE_ILUCreate(&solver_);
HYPRE_ILUSetType(solver_, get_optional(node, "ilu_type", 0));
HYPRE_ILUSetMaxIter(solver_, get_optional(node, "max_iterations", 20));
HYPRE_ILUSetTol(solver_, get_optional(node, "tolerance", 0.0));
HYPRE_ILUSetLocalReordering(solver_, get_optional(node,"local_reordering", 0));
HYPRE_ILUSetPrintLevel(solver_, get_optional(node, "print_level", 4));

// ILUK parameters
HYPRE_ILUSetLevelOfFill(solver_, get_optional(node, "fill", 0));

// ILUT parameters
HYPRE_ILUSetMaxNnzPerRow(solver_,get_optional(node, "max_nnz_per_row", 1000));
HYPRE_ILUSetDropThreshold(solver_,get_optional(node, "drop_threshold", 1.0e-2));

// 0 : Non-iterative algorithm (default)
// 1 : Asynchronous with in-place storage
// 2 : Asynchronous with explicit storage splitting
// 3 : Synchronous with explicit storage splitting
// 4 : Semi-synchronous with explicit storage splitting
// Iterative ILU is available only for zero fill-in and it depends on rocSparse
HYPRE_ILUSetIterativeSetupType(solver_,
get_optional(node, "algorithm_type", 0));
HYPRE_ILUSetIterativeSetupMaxIter(
solver_, get_optional(node, "max_ilu_iterations", 1));
HYPRE_ILUSetIterativeSetupTolerance(
solver_, get_optional(node, "iterative_ilu_tolerance", 1e-5));
// 0: iterative
// 1: direct (default)
HYPRE_ILUSetTriSolve(solver_, get_optional(node, "trisolve", 1));
// Jacobi iterations for lower and upper triangular solves
HYPRE_ILUSetLowerJacobiIters(solver_, get_optional(node, "lower_jacobi_iters", 5));
HYPRE_ILUSetUpperJacobiIters(solver_, get_optional(node, "upper_jacobi_iters", 5));

solverDestroyPtr_ = &HYPRE_ILUDestroy;
solverSetupPtr_ = &HYPRE_ILUSetup;
solverPrecondPtr_ = nullptr;
solverSolvePtr_ = &HYPRE_ILUSolve;
}

void HypreSystem::destroy_system() {
if (mat_)
HYPRE_IJMatrixDestroy(mat_);
Expand Down
5 changes: 5 additions & 0 deletions src/HypreSystem.h
Original file line number Diff line number Diff line change
Expand Up @@ -157,6 +157,11 @@ class HypreSystem {
void setup_bicg();
void setup_cg();

//! Setup ILU
void setup_ilu_precond();

void setup_ilu();

//! MPI Communicator object
MPI_Comm comm_;

Expand Down