An interface class to the suite of Hypre solvers and preconditioners to allow use of: More...
#include <hypre_solver.h>
Public Types | |
enum | Hypre_method_types { CG , GMRES , BiCGStab , BoomerAMG , Euclid , ParaSails , None } |
Enumerated flag to define which Hypre methods are used CAREFUL: DON'T CHANGE THE ORDER OF THESE! More... | |
Public Member Functions | |
HypreInterface () | |
Constructor. | |
~HypreInterface () | |
Destructor. | |
HypreInterface (const HypreInterface &)=delete | |
Broken copy constructor. | |
void | operator= (const HypreInterface &)=delete |
Broken assignment operator. | |
void | enable_hypre_error_messages () |
Turn on the hypre error messages. | |
void | disable_hypre_error_messages () |
Turn off hypre error messages. | |
unsigned | existing_solver () |
Function to return value of which solver (if any) is currently stored. | |
unsigned | existing_preconditioner () |
Function return value of which preconditioner (if any) is stored. | |
Protected Member Functions | |
void | hypre_clean_up_memory () |
Function deletes all solver data. | |
void | hypre_matrix_setup (CRDoubleMatrix *matrix_pt) |
Function which sets values of First_global_row, Last_global_row and other partitioning data and creates the distributed Hypre matrix (stored in Matrix_ij/Matrix_par) from the CRDoubleMatrix. | |
void | hypre_solver_setup () |
Sets up the data required for to use as an oomph-lib LinearSolver or Preconditioner. This must be called after the Hypre matrix has been generated using hypre_matrix_setup(...). | |
void | hypre_solve (const DoubleVector &rhs, DoubleVector &solution) |
Helper function performs a solve if any solver exists. | |
void | doc_hypre_parameters () |
Doc parameters used in hypre solver. | |
Protected Attributes | |
bool | Output_info |
Flag is true to output info and results of timings. | |
unsigned | Max_iter |
Maximum number of iterations used in solver. | |
double | Tolerance |
Tolerance used to terminate solver. | |
unsigned | Hypre_method |
Hypre method flag. Valid values are specified in enumeration. | |
unsigned | Internal_preconditioner |
Preconditioner method flag used with Hypre's PCG, GMRES and BiCGStab in solve(...) or resolve(...). Valid values are BoomerAMG, Euclid, ParaSails or None (all enumerated above), for any other value no preconditioner is set. | |
unsigned | AMG_print_level |
Used to set the Hypre printing level for AMG 0: no printout 1: print setup information 2: print solve information 3: print setup and solve information. | |
unsigned | AMG_max_levels |
Maximum number of levels used in AMG. | |
double | AMG_max_row_sum |
Parameter to identify diagonally dominant parts of the matrix in AMG. | |
bool | AMG_using_simple_smoothing |
Flag to determine whether simple smoothers (determined by the AMG_simple_smoother flag) or complex smoothers (determined by the AMG_complex_smoother flag are used in AMG. | |
unsigned | AMG_simple_smoother |
Simple smoothing methods used in BoomerAMG. To use these methods set AMG_using_simple_smoothing to true. Here are the current options (from hypre's HYPRE_parcsr_ls.h): | |
unsigned | AMG_complex_smoother |
Complex smoothing methods used in BoomerAMG. To use these methods set AMG_using_simple_smoothing to false Here are the current options (from hypre's HYPRE_parcsr_ls.h): | |
unsigned | AMG_smoother_iterations |
The number of smoother iterations to apply. | |
double | AMG_damping |
Damping factor for BoomerAMG smoothed Jacobi or hybrid SOR. | |
double | AMG_strength |
Connection strength threshold parameter for BoomerAMG. | |
double | AMG_truncation |
Interpolation truncation factor for BoomerAMG. | |
unsigned | AMG_coarsening |
AMG coarsening strategy. Coarsening types include: 0 = CLJP (parallel coarsening using independent sets) 1 = classical RS with no boundary treatment (not recommended in parallel) 3 = modified RS with 3rd pass to add C points on the boundaries 6 = Falgout (uses 1 then CLJP using interior coarse points as first independent set) 8 = PMIS (parallel coarsening using independent sets - lower complexities than 0, maybe also slower convergence) 10= HMIS (one pass RS on each processor then PMIS on interior coarse points as first independent set) 11= One pass RS on each processor (not recommended) | |
int | ParaSails_symmetry |
ParaSails symmetry flag, used to inform ParaSails of Symmetry of definitenss of problem and type of ParaSails preconditioner: 0 = nonsymmetric and/or indefinite problem, nonsymmetric preconditioner 1 = SPD problem, and SPD (factored preconditioner) 2 = nonsymmetric, definite problem and SDP (factored preconditoner) | |
int | ParaSails_nlevel |
ParaSails nlevel parameter. | |
double | ParaSails_thresh |
ParaSails thresh parameter. | |
double | ParaSails_filter |
ParaSails filter parameter. | |
double | Euclid_droptol |
Euclid drop tolerance for ILU(k) and ILUT factorization. | |
bool | Euclid_rowScale |
Flag to switch on Euclid row scaling. | |
bool | Euclid_using_ILUT |
Flag to determine if ILUT (if true) or ILU(k) is used in Euclid. | |
bool | Euclid_using_BJ |
Flag to determine if Block Jacobi is used instead of PILU. | |
int | Euclid_level |
Euclid level parameter for ILU(k) factorization. | |
unsigned | Euclid_print_level |
Flag to set the level of printing from Euclid when the Euclid destructor is called 0: no printing (default) 1: prints summary of runtime settings and timings 2: as 1 plus prints memory usage. | |
bool | AMGEuclidSmoother_use_block_jacobi |
bool | AMGEuclidSmoother_use_row_scaling |
bool | AMGEuclidSmoother_use_ilut |
unsigned | AMGEuclidSmoother_level |
double | AMGEuclidSmoother_drop_tol |
unsigned | AMGEuclidSmoother_print_level |
unsigned | Krylov_print_level |
Used to set the Hypre printing level for the Krylov subspace solvers. | |
bool | Hypre_error_messages |
Flag to determine if non-zero values of the Hypre error flag plus Hypre error messages are output to screen at various points in the solve function, i.e. after: | |
bool | Delete_input_data |
Internal flag which is true when hypre_setup or hypre_solve can delete input matrix. | |
bool | Using_distributed_rhs |
Internal flag which tell the solver if the rhs Vector is distributed or not. | |
bool | Returning_distributed_solution |
Internal flag which tell the solver if the solution Vector to be returned is distributed or not. | |
Private Attributes | |
bool | Delete_matrix |
Hypre copies matrix data from oomph-lib's CRDoubleMatrix into its own data structures, doubling the memory requirements. As far as the Hypre solvers are concerned the oomph-lib matrix is no longer required and could be deleted to save memory. However, this will not be the expected behaviour and therefore needs to be requested explicitly by the user by changing this flag from false (its default) to true. | |
HYPRE_IJMatrix | Matrix_ij |
The Hypre_IJMatrix version of the matrix used in solve(...), resolve(...) or preconditioner_solve(...). | |
HYPRE_ParCSRMatrix | Matrix_par |
The Hypre_ParCSRMatrix version of the matrix used in solve(...), resolve(...) or preconditioner_solve(...). | |
HYPRE_Solver | Solver |
The Hypre solver used in solve(...), resolve(...) or preconditioner_solve(...). [This is a C structure!]. | |
HYPRE_Solver | Preconditioner |
The internal Hypre preconditioner used in conjunction with Solver. [This is a C structure!]. | |
unsigned | Existing_solver |
Used to keep track of which solver (if any) is currently stored. | |
unsigned | Existing_preconditioner |
Used to keep track of which preconditioner (if any) is currently stored. | |
LinearAlgebraDistribution * | Hypre_distribution_pt |
the distribution for this helpers- | |
An interface class to the suite of Hypre solvers and preconditioners to allow use of:
BoomerAMG (AMG), CG, GMRES or BiCGStab, Euclid (ILU) or ParaSails (Approximate inverse)
Hypre's Krylov subspace solvers (CG, GMRES and BiCGStab) may be preconditioned using:
BoomerAMG, Euclid or ParaSails
Definition at line 138 of file hypre_solver.h.
Enumerated flag to define which Hypre methods are used CAREFUL: DON'T CHANGE THE ORDER OF THESE!
Enumerator | |
---|---|
CG | |
GMRES | |
BiCGStab | |
BoomerAMG | |
Euclid | |
ParaSails | |
None |
Definition at line 286 of file hypre_solver.h.
|
inline |
Constructor.
Definition at line 142 of file hypre_solver.h.
References oomph::HypreHelpers::AMG_coarsening, AMG_coarsening, AMG_complex_smoother, AMG_damping, AMG_max_levels, AMG_max_row_sum, AMG_print_level, AMG_simple_smoother, AMG_smoother_iterations, oomph::HypreHelpers::AMG_strength, AMG_strength, oomph::HypreHelpers::AMG_truncation, AMG_truncation, AMG_using_simple_smoothing, AMGEuclidSmoother_drop_tol, AMGEuclidSmoother_level, AMGEuclidSmoother_print_level, AMGEuclidSmoother_use_block_jacobi, AMGEuclidSmoother_use_ilut, AMGEuclidSmoother_use_row_scaling, oomph::MPI_Helpers::communicator_pt(), e, Euclid_droptol, Euclid_level, Euclid_print_level, Euclid_rowScale, Euclid_using_BJ, Euclid_using_ILUT, Existing_preconditioner, Existing_solver, GMRES, Hypre_distribution_pt, Hypre_error_messages, Hypre_method, Internal_preconditioner, Krylov_print_level, Max_iter, oomph::MPI_Helpers::mpi_has_been_initialised(), None, oomph::HypreHelpers::Number_of_active_hypre_solvers, Output_info, ParaSails_filter, ParaSails_nlevel, ParaSails_symmetry, ParaSails_thresh, and Tolerance.
|
inline |
Destructor.
Definition at line 238 of file hypre_solver.h.
References hypre_clean_up_memory(), Hypre_distribution_pt, and oomph::HypreHelpers::Number_of_active_hypre_solvers.
|
delete |
Broken copy constructor.
|
inline |
Turn off hypre error messages.
Definition at line 279 of file hypre_solver.h.
References Hypre_error_messages.
|
inlineprotected |
Doc parameters used in hypre solver.
Definition at line 514 of file hypre_solver.h.
References AMG_coarsening, AMG_complex_smoother, AMG_damping, AMG_max_levels, AMG_max_row_sum, AMG_print_level, AMG_simple_smoother, AMG_smoother_iterations, AMG_strength, AMG_truncation, AMG_using_simple_smoothing, Hypre_method, Internal_preconditioner, Krylov_print_level, Max_iter, oomph::oomph_info, and Tolerance.
|
inline |
Turn on the hypre error messages.
Definition at line 273 of file hypre_solver.h.
References Hypre_error_messages.
|
inline |
Function return value of which preconditioner (if any) is stored.
Definition at line 304 of file hypre_solver.h.
References Existing_preconditioner.
|
inline |
Function to return value of which solver (if any) is currently stored.
Definition at line 298 of file hypre_solver.h.
References Existing_solver.
Referenced by oomph::HyprePreconditioner::preconditioner_solve(), and oomph::HypreSolver::resolve().
|
protected |
Function deletes all solver data.
hypre_clean_up_memory() deletes any existing Hypre solver and Hypre matrix
Definition at line 1364 of file hypre_solver.cc.
References BoomerAMG, oomph::HypreHelpers::check_HYPRE_error_flag(), Euclid, Existing_preconditioner, Existing_solver, Hypre_error_messages, Matrix_ij, None, ParaSails, and Solver.
Referenced by oomph::HypreSolver::clean_up_memory(), oomph::HyprePreconditioner::clean_up_memory(), and ~HypreInterface().
|
protected |
Function which sets values of First_global_row, Last_global_row and other partitioning data and creates the distributed Hypre matrix (stored in Matrix_ij/Matrix_par) from the CRDoubleMatrix.
Helper function which creates a Hypre matrix from a CRDoubleMatrix If OOMPH-LIB has been set up for MPI use, the Hypre matrix is distributed over the available processors.
Definition at line 548 of file hypre_solver.cc.
References oomph::HypreHelpers::check_HYPRE_error_flag(), oomph::CRDoubleMatrix::clear(), oomph::LinearAlgebraDistribution::communicator_pt(), oomph::HypreHelpers::create_HYPRE_Matrix(), Delete_input_data, oomph::DistributableLinearAlgebraObject::distribution_pt(), hypre__global_error, Hypre_distribution_pt, Hypre_error_messages, Matrix_ij, Matrix_par, oomph::CRDoubleMatrix::nrow(), and oomph::oomph_info.
Referenced by oomph::HyprePreconditioner::setup(), oomph::HypreSolver::solve(), and oomph::HypreSolver::solve().
|
protected |
Helper function performs a solve if any solver exists.
Helper function performs a solve if solver data has been set up using hypre_solver_setup(...).
Definition at line 1149 of file hypre_solver.cc.
References BoomerAMG, oomph::HypreHelpers::check_HYPRE_error_flag(), oomph::MPI_Helpers::communicator_pt(), oomph::HypreHelpers::create_HYPRE_Vector(), Euclid, Existing_solver, oomph::LinearAlgebraDistribution::first_row(), Hypre_distribution_pt, Hypre_error_messages, Hypre_method, i, Matrix_par, oomph::MPI_Helpers::mpi_has_been_initialised(), oomph::LinearAlgebraDistribution::nrow_local(), oomph::oomph_info, Output_info, ParaSails, Solver, and oomph::TimingHelpers::timer().
Referenced by oomph::HyprePreconditioner::preconditioner_solve(), oomph::HypreSolver::resolve(), oomph::HypreSolver::solve(), and oomph::HypreSolver::solve().
|
protected |
Sets up the data required for to use as an oomph-lib LinearSolver or Preconditioner. This must be called after the Hypre matrix has been generated using hypre_matrix_setup(...).
Sets up the solver data required for use in an oomph-lib LinearSolver or Preconditioner, once the Hypre matrix has been generated using hypre_matrix_setup(...).
Definition at line 592 of file hypre_solver.cc.
References AMG_coarsening, AMG_complex_smoother, AMG_damping, AMG_max_levels, AMG_max_row_sum, AMG_print_level, AMG_simple_smoother, AMG_smoother_iterations, AMG_strength, AMG_truncation, AMG_using_simple_smoothing, AMGEuclidSmoother_drop_tol, AMGEuclidSmoother_level, AMGEuclidSmoother_print_level, AMGEuclidSmoother_use_block_jacobi, AMGEuclidSmoother_use_ilut, AMGEuclidSmoother_use_row_scaling, BiCGStab, BoomerAMG, CG, oomph::HypreHelpers::check_HYPRE_error_flag(), oomph::LinearAlgebraDistribution::communicator_pt(), oomph::HypreHelpers::create_HYPRE_Vector(), Euclid, Euclid_droptol, Euclid_level, Euclid_print_level, Euclid_rowScale, oomph::HypreHelpers::euclid_settings_helper(), Euclid_using_BJ, Euclid_using_ILUT, Existing_preconditioner, Existing_solver, GMRES, hypre__global_error, Hypre_distribution_pt, Hypre_error_messages, Hypre_method, Internal_preconditioner, Krylov_print_level, Matrix_par, Max_iter, oomph::oomph_info, Output_info, ParaSails, ParaSails_filter, ParaSails_nlevel, ParaSails_symmetry, ParaSails_thresh, Solver, oomph::TimingHelpers::timer(), and Tolerance.
Referenced by oomph::HyprePreconditioner::setup(), oomph::HypreSolver::solve(), and oomph::HypreSolver::solve().
|
delete |
Broken assignment operator.
|
protected |
AMG coarsening strategy. Coarsening types include: 0 = CLJP (parallel coarsening using independent sets) 1 = classical RS with no boundary treatment (not recommended in parallel) 3 = modified RS with 3rd pass to add C points on the boundaries 6 = Falgout (uses 1 then CLJP using interior coarse points as first independent set) 8 = PMIS (parallel coarsening using independent sets - lower complexities than 0, maybe also slower convergence) 10= HMIS (one pass RS on each processor then PMIS on interior coarse points as first independent set) 11= One pass RS on each processor (not recommended)
Definition at line 459 of file hypre_solver.h.
Referenced by oomph::HypreSolver::amg_coarsening(), oomph::HyprePreconditioner::amg_coarsening(), doc_hypre_parameters(), hypre_solver_setup(), and HypreInterface().
|
protected |
Complex smoothing methods used in BoomerAMG. To use these methods set AMG_using_simple_smoothing to false Here are the current options (from hypre's HYPRE_parcsr_ls.h):
Definition at line 433 of file hypre_solver.h.
Referenced by oomph::HypreSolver::amg_complex_smoother(), oomph::HyprePreconditioner::amg_complex_smoother(), doc_hypre_parameters(), hypre_solver_setup(), and HypreInterface().
|
protected |
Damping factor for BoomerAMG smoothed Jacobi or hybrid SOR.
Definition at line 439 of file hypre_solver.h.
Referenced by oomph::HypreSolver::amg_damping(), oomph::HyprePreconditioner::amg_damping(), doc_hypre_parameters(), hypre_solver_setup(), and HypreInterface().
|
protected |
Maximum number of levels used in AMG.
Definition at line 354 of file hypre_solver.h.
Referenced by oomph::HypreSolver::amg_max_levels(), oomph::HyprePreconditioner::amg_max_levels(), doc_hypre_parameters(), hypre_solver_setup(), and HypreInterface().
|
protected |
Parameter to identify diagonally dominant parts of the matrix in AMG.
Definition at line 357 of file hypre_solver.h.
Referenced by oomph::HypreSolver::amg_max_row_sum(), oomph::HyprePreconditioner::amg_max_row_sum(), doc_hypre_parameters(), hypre_solver_setup(), and HypreInterface().
|
protected |
Used to set the Hypre printing level for AMG 0: no printout 1: print setup information 2: print solve information 3: print setup and solve information.
Definition at line 351 of file hypre_solver.h.
Referenced by oomph::HypreSolver::amg_print_level(), oomph::HyprePreconditioner::amg_print_level(), doc_hypre_parameters(), hypre_solver_setup(), and HypreInterface().
|
protected |
Simple smoothing methods used in BoomerAMG. To use these methods set AMG_using_simple_smoothing to true. Here are the current options (from hypre's HYPRE_parcsr_ls.h):
(Optional) Defines the smoother to be used. It uses the given smoother on the fine grid, the up and the down cycle and sets the solver on the coarsest level to Gaussian elimination (9). The default is -Gauss-Seidel, forward solve (13) on the down cycle and backward solve (14) on the up cycle.
There are the following options for relax_type:
Definition at line 411 of file hypre_solver.h.
Referenced by oomph::HypreSolver::amg_simple_smoother(), oomph::HyprePreconditioner::amg_simple_smoother(), doc_hypre_parameters(), hypre_solver_setup(), and HypreInterface().
|
protected |
The number of smoother iterations to apply.
Definition at line 436 of file hypre_solver.h.
Referenced by oomph::HypreSolver::amg_smoother_iterations(), oomph::HyprePreconditioner::amg_smoother_iterations(), doc_hypre_parameters(), hypre_solver_setup(), and HypreInterface().
|
protected |
Connection strength threshold parameter for BoomerAMG.
Definition at line 442 of file hypre_solver.h.
Referenced by oomph::HypreSolver::amg_strength(), oomph::HyprePreconditioner::amg_strength(), doc_hypre_parameters(), hypre_solver_setup(), and HypreInterface().
|
protected |
Interpolation truncation factor for BoomerAMG.
Definition at line 445 of file hypre_solver.h.
Referenced by oomph::HypreSolver::amg_truncation(), oomph::HyprePreconditioner::amg_truncation(), doc_hypre_parameters(), hypre_solver_setup(), and HypreInterface().
|
protected |
Flag to determine whether simple smoothers (determined by the AMG_simple_smoother flag) or complex smoothers (determined by the AMG_complex_smoother flag are used in AMG.
Definition at line 362 of file hypre_solver.h.
Referenced by oomph::HypreSolver::amg_using_complex_smoothing(), oomph::HyprePreconditioner::amg_using_complex_smoothing(), oomph::HypreSolver::amg_using_simple_smoothing(), oomph::HyprePreconditioner::amg_using_simple_smoothing(), oomph::HyprePreconditioner::amg_using_simple_smoothing_flag(), doc_hypre_parameters(), hypre_solver_setup(), and HypreInterface().
|
protected |
Definition at line 506 of file hypre_solver.h.
Referenced by hypre_solver_setup(), and HypreInterface().
|
protected |
Definition at line 505 of file hypre_solver.h.
Referenced by hypre_solver_setup(), and HypreInterface().
|
protected |
Definition at line 507 of file hypre_solver.h.
Referenced by hypre_solver_setup(), and HypreInterface().
|
protected |
Definition at line 502 of file hypre_solver.h.
Referenced by hypre_solver_setup(), and HypreInterface().
|
protected |
Definition at line 504 of file hypre_solver.h.
Referenced by hypre_solver_setup(), and HypreInterface().
|
protected |
Definition at line 503 of file hypre_solver.h.
Referenced by hypre_solver_setup(), and HypreInterface().
|
protected |
Internal flag which is true when hypre_setup or hypre_solve can delete input matrix.
Definition at line 634 of file hypre_solver.h.
Referenced by hypre_matrix_setup(), oomph::HyprePreconditioner::setup(), oomph::HypreSolver::solve(), and oomph::HypreSolver::solve().
|
private |
Hypre copies matrix data from oomph-lib's CRDoubleMatrix into its own data structures, doubling the memory requirements. As far as the Hypre solvers are concerned the oomph-lib matrix is no longer required and could be deleted to save memory. However, this will not be the expected behaviour and therefore needs to be requested explicitly by the user by changing this flag from false (its default) to true.
Definition at line 654 of file hypre_solver.h.
|
protected |
Euclid drop tolerance for ILU(k) and ILUT factorization.
Definition at line 479 of file hypre_solver.h.
Referenced by oomph::HypreSolver::euclid_droptol(), oomph::HyprePreconditioner::euclid_droptol(), hypre_solver_setup(), and HypreInterface().
|
protected |
Euclid level parameter for ILU(k) factorization.
Definition at line 491 of file hypre_solver.h.
Referenced by oomph::HypreSolver::euclid_level(), oomph::HyprePreconditioner::euclid_level(), hypre_solver_setup(), and HypreInterface().
|
protected |
Flag to set the level of printing from Euclid when the Euclid destructor is called 0: no printing (default) 1: prints summary of runtime settings and timings 2: as 1 plus prints memory usage.
Definition at line 498 of file hypre_solver.h.
Referenced by oomph::HypreSolver::euclid_print_level(), oomph::HyprePreconditioner::euclid_print_level(), hypre_solver_setup(), and HypreInterface().
|
protected |
Flag to switch on Euclid row scaling.
Definition at line 482 of file hypre_solver.h.
Referenced by oomph::HypreSolver::disable_euclid_rowScale(), oomph::HyprePreconditioner::disable_euclid_rowScale(), oomph::HypreSolver::enable_euclid_rowScale(), oomph::HyprePreconditioner::enable_euclid_rowScale(), hypre_solver_setup(), and HypreInterface().
|
protected |
Flag to determine if Block Jacobi is used instead of PILU.
Definition at line 488 of file hypre_solver.h.
Referenced by oomph::HypreSolver::disable_euclid_using_BJ(), oomph::HyprePreconditioner::disable_euclid_using_BJ(), oomph::HypreSolver::enable_euclid_using_BJ(), oomph::HyprePreconditioner::enable_euclid_using_BJ(), hypre_solver_setup(), and HypreInterface().
|
protected |
Flag to determine if ILUT (if true) or ILU(k) is used in Euclid.
Definition at line 485 of file hypre_solver.h.
Referenced by oomph::HypreSolver::euclid_using_ILUK(), oomph::HyprePreconditioner::euclid_using_ILUK(), oomph::HypreSolver::euclid_using_ILUT(), oomph::HyprePreconditioner::euclid_using_ILUT(), hypre_solver_setup(), and HypreInterface().
|
private |
Used to keep track of which preconditioner (if any) is currently stored.
Definition at line 676 of file hypre_solver.h.
Referenced by existing_preconditioner(), hypre_clean_up_memory(), hypre_solver_setup(), and HypreInterface().
|
private |
Used to keep track of which solver (if any) is currently stored.
Definition at line 673 of file hypre_solver.h.
Referenced by existing_solver(), hypre_clean_up_memory(), hypre_solve(), hypre_solver_setup(), and HypreInterface().
|
private |
the distribution for this helpers-
Definition at line 679 of file hypre_solver.h.
Referenced by hypre_matrix_setup(), hypre_solve(), hypre_solver_setup(), HypreInterface(), and ~HypreInterface().
|
protected |
Flag to determine if non-zero values of the Hypre error flag plus Hypre error messages are output to screen at various points in the solve function, i.e. after:
Definition at line 630 of file hypre_solver.h.
Referenced by disable_hypre_error_messages(), enable_hypre_error_messages(), hypre_clean_up_memory(), hypre_matrix_setup(), hypre_solve(), hypre_solver_setup(), and HypreInterface().
|
protected |
Hypre method flag. Valid values are specified in enumeration.
Definition at line 338 of file hypre_solver.h.
Referenced by doc_hypre_parameters(), oomph::HypreSolver::hypre_method(), oomph::HyprePreconditioner::hypre_method(), hypre_solve(), hypre_solver_setup(), HypreInterface(), oomph::HyprePreconditioner::HyprePreconditioner(), oomph::HyprePreconditioner::use_BoomerAMG(), oomph::HyprePreconditioner::use_Euclid(), and oomph::HyprePreconditioner::use_ParaSails().
|
protected |
Preconditioner method flag used with Hypre's PCG, GMRES and BiCGStab in solve(...) or resolve(...). Valid values are BoomerAMG, Euclid, ParaSails or None (all enumerated above), for any other value no preconditioner is set.
Definition at line 344 of file hypre_solver.h.
Referenced by doc_hypre_parameters(), hypre_solver_setup(), HypreInterface(), oomph::HyprePreconditioner::HyprePreconditioner(), oomph::HypreSolver::internal_preconditioner(), and oomph::HyprePreconditioner::internal_preconditioner().
|
protected |
Used to set the Hypre printing level for the Krylov subspace solvers.
Definition at line 511 of file hypre_solver.h.
Referenced by doc_hypre_parameters(), hypre_solver_setup(), HypreInterface(), and oomph::HypreSolver::krylov_print_level().
|
private |
The Hypre_IJMatrix version of the matrix used in solve(...), resolve(...) or preconditioner_solve(...).
Definition at line 658 of file hypre_solver.h.
Referenced by hypre_clean_up_memory(), and hypre_matrix_setup().
|
private |
The Hypre_ParCSRMatrix version of the matrix used in solve(...), resolve(...) or preconditioner_solve(...).
Definition at line 662 of file hypre_solver.h.
Referenced by hypre_matrix_setup(), hypre_solve(), and hypre_solver_setup().
|
protected |
Maximum number of iterations used in solver.
Definition at line 332 of file hypre_solver.h.
Referenced by oomph::HyprePreconditioner::amg_iterations(), doc_hypre_parameters(), hypre_solver_setup(), HypreInterface(), oomph::HyprePreconditioner::HyprePreconditioner(), oomph::HypreSolver::max_iter(), and oomph::HyprePreconditioner::set_amg_iterations().
|
protected |
Flag is true to output info and results of timings.
Definition at line 329 of file hypre_solver.h.
Referenced by hypre_solve(), hypre_solver_setup(), HypreInterface(), oomph::HyprePreconditioner::preconditioner_solve(), oomph::HypreSolver::resolve(), oomph::HyprePreconditioner::setup(), oomph::HypreSolver::solve(), and oomph::HypreSolver::solve().
|
protected |
ParaSails filter parameter.
Definition at line 476 of file hypre_solver.h.
Referenced by hypre_solver_setup(), HypreInterface(), oomph::HypreSolver::parasails_filter(), and oomph::HyprePreconditioner::parasails_filter().
|
protected |
ParaSails nlevel parameter.
Definition at line 470 of file hypre_solver.h.
Referenced by hypre_solver_setup(), HypreInterface(), oomph::HypreSolver::parasails_nlevel(), and oomph::HyprePreconditioner::parasails_nlevel().
|
protected |
ParaSails symmetry flag, used to inform ParaSails of Symmetry of definitenss of problem and type of ParaSails preconditioner: 0 = nonsymmetric and/or indefinite problem, nonsymmetric preconditioner 1 = SPD problem, and SPD (factored preconditioner) 2 = nonsymmetric, definite problem and SDP (factored preconditoner)
Definition at line 467 of file hypre_solver.h.
Referenced by hypre_solver_setup(), HypreInterface(), oomph::HypreSolver::parasails_symmetry(), and oomph::HyprePreconditioner::parasails_symmetry().
|
protected |
ParaSails thresh parameter.
Definition at line 473 of file hypre_solver.h.
Referenced by hypre_solver_setup(), HypreInterface(), oomph::HypreSolver::parasails_thresh(), and oomph::HyprePreconditioner::parasails_thresh().
|
private |
The internal Hypre preconditioner used in conjunction with Solver. [This is a C structure!].
Definition at line 670 of file hypre_solver.h.
|
protected |
Internal flag which tell the solver if the solution Vector to be returned is distributed or not.
Definition at line 643 of file hypre_solver.h.
|
private |
The Hypre solver used in solve(...), resolve(...) or preconditioner_solve(...). [This is a C structure!].
Definition at line 666 of file hypre_solver.h.
Referenced by hypre_clean_up_memory(), hypre_solve(), and hypre_solver_setup().
|
protected |
Tolerance used to terminate solver.
Definition at line 335 of file hypre_solver.h.
Referenced by doc_hypre_parameters(), hypre_solver_setup(), HypreInterface(), oomph::HyprePreconditioner::HyprePreconditioner(), and oomph::HypreSolver::tolerance().
|
protected |
Internal flag which tell the solver if the rhs Vector is distributed or not.
Definition at line 639 of file hypre_solver.h.