Public Types | Public Member Functions | Public Attributes | Protected Member Functions | Protected Attributes | Private Attributes | List of all members
oomph::HypreInterface Class Reference

An interface class to the suite of Hypre solvers and preconditioners to allow use of: More...

#include <hypre_solver.h>

+ Inheritance diagram for oomph::HypreInterface:

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. More...
 
 ~HypreInterface ()
 Destructor. More...
 
 HypreInterface (const HypreInterface &)=delete
 Broken copy constructor. More...
 
void operator= (const HypreInterface &)=delete
 Broken assignment operator. More...
 
void enable_hypre_error_messages ()
 Turn on the hypre error messages. More...
 
void disable_hypre_error_messages ()
 Turn off hypre error messages. More...
 
unsigned existing_solver ()
 Function to return value of which solver (if any) is currently stored. More...
 
unsigned existing_preconditioner ()
 Function return value of which preconditioner (if any) is stored. More...
 

Public Attributes

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
 

Protected Member Functions

void hypre_clean_up_memory ()
 Function deletes all solver data. More...
 
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. More...
 
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(...). More...
 
void hypre_solve (const DoubleVector &rhs, DoubleVector &solution)
 Helper function performs a solve if any solver exists. More...
 

Protected Attributes

bool Output_info
 Flag is true to output info and results of timings. More...
 
unsigned Max_iter
 Maximum number of iterations used in solver. More...
 
double Tolerance
 Tolerance used to terminate solver. More...
 
unsigned Hypre_method
 Hypre method flag. Valid values are specified in enumeration. More...
 
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. More...
 
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. More...
 
unsigned AMG_max_levels
 Maximum number of levels used in AMG. More...
 
double AMG_max_row_sum
 Parameter to identify diagonally dominant parts of the matrix in AMG. More...
 
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. More...
 
unsigned AMG_simple_smoother
 Simple smoothing methods used in BoomerAMG. Relaxation types include: 0 = Jacobi 1 = Gauss-Seidel, sequential (very slow in parallel!) 2 = Gauss-Seidel, interior points in parallel, boundary sequential (slow in parallel!) 3 = hybrid Gauss-Seidel or SOR, forward solve 4 = hybrid Gauss-Seidel or SOR, backward solve 6 = hybrid symmetric Gauss-Seidel or SSOR To use these methods set AMG_using_simple_smoothing to true. More...
 
unsigned AMG_complex_smoother
 Complex smoothing methods used in BoomerAMG. Relaxation types are: 6 = Schwarz 7 = Pilut 8 = ParaSails 9 = Euclid To use these methods set AMG_using_simple_smoothing to false. More...
 
unsigned AMG_smoother_iterations
 The number of smoother iterations to apply. More...
 
double AMG_damping
 Damping factor for BoomerAMG smoothed Jacobi or hybrid SOR. More...
 
double AMG_strength
 Connection strength threshold parameter for BoomerAMG. More...
 
double AMG_truncation
 Interpolation truncation factor for BoomerAMG. More...
 
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) More...
 
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) More...
 
int ParaSails_nlevel
 ParaSails nlevel parameter. More...
 
double ParaSails_thresh
 ParaSails thresh parameter. More...
 
double ParaSails_filter
 ParaSails filter parameter. More...
 
double Euclid_droptol
 Euclid drop tolerance for ILU(k) and ILUT factorization. More...
 
bool Euclid_rowScale
 Flag to switch on Euclid row scaling. More...
 
bool Euclid_using_ILUT
 Flag to determine if ILUT (if true) or ILU(k) is used in Euclid. More...
 
bool Euclid_using_BJ
 Flag to determine if Block Jacobi is used instead of PILU. More...
 
int Euclid_level
 Euclid level parameter for ILU(k) factorization. More...
 
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. More...
 
unsigned Krylov_print_level
 Used to set the Hypre printing level for the Krylov subspace solvers. More...
 
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: More...
 
bool Delete_input_data
 Internal flag which is true when hypre_setup or hypre_solve can delete input matrix. More...
 
bool Using_distributed_rhs
 Internal flag which tell the solver if the rhs Vector is distributed or not. More...
 
bool Returning_distributed_solution
 Internal flag which tell the solver if the solution Vector to be returned is distributed or not. More...
 

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. More...
 
HYPRE_IJMatrix Matrix_ij
 The Hypre_IJMatrix version of the matrix used in solve(...), resolve(...) or preconditioner_solve(...). More...
 
HYPRE_ParCSRMatrix Matrix_par
 The Hypre_ParCSRMatrix version of the matrix used in solve(...), resolve(...) or preconditioner_solve(...). More...
 
HYPRE_Solver Solver
 The Hypre solver used in solve(...), resolve(...) or preconditioner_solve(...). [This is a C structure!]. More...
 
HYPRE_Solver Preconditioner
 The internal Hypre preconditioner used in conjunction with Solver. [This is a C structure!]. More...
 
unsigned Existing_solver
 Used to keep track of which solver (if any) is currently stored. More...
 
unsigned Existing_preconditioner
 Used to keep track of which preconditioner (if any) is currently stored. More...
 
LinearAlgebraDistributionHypre_distribution_pt
 the distribution for this helpers- More...
 

Detailed Description

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 133 of file hypre_solver.h.

Member Enumeration Documentation

◆ Hypre_method_types

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 255 of file hypre_solver.h.

Constructor & Destructor Documentation

◆ HypreInterface() [1/2]

oomph::HypreInterface::HypreInterface ( )
inline

◆ ~HypreInterface()

oomph::HypreInterface::~HypreInterface ( )
inline

Destructor.

Definition at line 226 of file hypre_solver.h.

References hypre_clean_up_memory(), and Hypre_distribution_pt.

◆ HypreInterface() [2/2]

oomph::HypreInterface::HypreInterface ( const HypreInterface )
delete

Broken copy constructor.

Member Function Documentation

◆ disable_hypre_error_messages()

void oomph::HypreInterface::disable_hypre_error_messages ( )
inline

Turn off hypre error messages.

Definition at line 248 of file hypre_solver.h.

References Hypre_error_messages.

◆ enable_hypre_error_messages()

void oomph::HypreInterface::enable_hypre_error_messages ( )
inline

Turn on the hypre error messages.

Definition at line 242 of file hypre_solver.h.

References Hypre_error_messages.

◆ existing_preconditioner()

unsigned oomph::HypreInterface::existing_preconditioner ( )
inline

Function return value of which preconditioner (if any) is stored.

Definition at line 273 of file hypre_solver.h.

References Existing_preconditioner.

◆ existing_solver()

unsigned oomph::HypreInterface::existing_solver ( )
inline

Function to return value of which solver (if any) is currently stored.

Definition at line 267 of file hypre_solver.h.

References Existing_solver.

Referenced by oomph::HyprePreconditioner::preconditioner_solve(), and oomph::HypreSolver::resolve().

◆ hypre_clean_up_memory()

void oomph::HypreInterface::hypre_clean_up_memory ( )
protected

◆ hypre_matrix_setup()

void oomph::HypreInterface::hypre_matrix_setup ( CRDoubleMatrix matrix_pt)
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 513 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(), and oomph::HypreSolver::solve().

◆ hypre_solve()

void oomph::HypreInterface::hypre_solve ( const DoubleVector rhs,
DoubleVector solution 
)
protected

◆ hypre_solver_setup()

void oomph::HypreInterface::hypre_solver_setup ( )
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 558 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, i, 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(), and oomph::HypreSolver::solve().

◆ operator=()

void oomph::HypreInterface::operator= ( const HypreInterface )
delete

Broken assignment operator.

Member Data Documentation

◆ AMG_coarsening

unsigned oomph::HypreInterface::AMG_coarsening
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 386 of file hypre_solver.h.

Referenced by oomph::HypreSolver::amg_coarsening(), oomph::HyprePreconditioner::amg_coarsening(), hypre_solver_setup(), and HypreInterface().

◆ AMG_complex_smoother

unsigned oomph::HypreInterface::AMG_complex_smoother
protected

Complex smoothing methods used in BoomerAMG. Relaxation types are: 6 = Schwarz 7 = Pilut 8 = ParaSails 9 = Euclid To use these methods set AMG_using_simple_smoothing to false.

Definition at line 360 of file hypre_solver.h.

Referenced by oomph::HypreSolver::amg_complex_smoother(), oomph::HyprePreconditioner::amg_complex_smoother(), hypre_solver_setup(), and HypreInterface().

◆ AMG_damping

double oomph::HypreInterface::AMG_damping
protected

Damping factor for BoomerAMG smoothed Jacobi or hybrid SOR.

Definition at line 366 of file hypre_solver.h.

Referenced by oomph::HypreSolver::amg_damping(), oomph::HyprePreconditioner::amg_damping(), hypre_solver_setup(), and HypreInterface().

◆ AMG_max_levels

unsigned oomph::HypreInterface::AMG_max_levels
protected

Maximum number of levels used in AMG.

Definition at line 330 of file hypre_solver.h.

Referenced by oomph::HypreSolver::amg_max_levels(), oomph::HyprePreconditioner::amg_max_levels(), hypre_solver_setup(), and HypreInterface().

◆ AMG_max_row_sum

double oomph::HypreInterface::AMG_max_row_sum
protected

Parameter to identify diagonally dominant parts of the matrix in AMG.

Definition at line 333 of file hypre_solver.h.

Referenced by oomph::HypreSolver::amg_max_row_sum(), oomph::HyprePreconditioner::amg_max_row_sum(), hypre_solver_setup(), and HypreInterface().

◆ AMG_print_level

unsigned oomph::HypreInterface::AMG_print_level
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 327 of file hypre_solver.h.

Referenced by oomph::HypreSolver::amg_print_level(), oomph::HyprePreconditioner::amg_print_level(), hypre_solver_setup(), and HypreInterface().

◆ AMG_simple_smoother

unsigned oomph::HypreInterface::AMG_simple_smoother
protected

Simple smoothing methods used in BoomerAMG. Relaxation types include: 0 = Jacobi 1 = Gauss-Seidel, sequential (very slow in parallel!) 2 = Gauss-Seidel, interior points in parallel, boundary sequential (slow in parallel!) 3 = hybrid Gauss-Seidel or SOR, forward solve 4 = hybrid Gauss-Seidel or SOR, backward solve 6 = hybrid symmetric Gauss-Seidel or SSOR To use these methods set AMG_using_simple_smoothing to true.

Definition at line 351 of file hypre_solver.h.

Referenced by oomph::HypreSolver::amg_simple_smoother(), oomph::HyprePreconditioner::amg_simple_smoother(), hypre_solver_setup(), and HypreInterface().

◆ AMG_smoother_iterations

unsigned oomph::HypreInterface::AMG_smoother_iterations
protected

The number of smoother iterations to apply.

Definition at line 363 of file hypre_solver.h.

Referenced by oomph::HypreSolver::amg_smoother_iterations(), oomph::HyprePreconditioner::amg_smoother_iterations(), hypre_solver_setup(), and HypreInterface().

◆ AMG_strength

double oomph::HypreInterface::AMG_strength
protected

Connection strength threshold parameter for BoomerAMG.

Definition at line 369 of file hypre_solver.h.

Referenced by oomph::HypreSolver::amg_strength(), oomph::HyprePreconditioner::amg_strength(), hypre_solver_setup(), and HypreInterface().

◆ AMG_truncation

double oomph::HypreInterface::AMG_truncation
protected

Interpolation truncation factor for BoomerAMG.

Definition at line 372 of file hypre_solver.h.

Referenced by oomph::HypreSolver::amg_truncation(), oomph::HyprePreconditioner::amg_truncation(), hypre_solver_setup(), and HypreInterface().

◆ AMG_using_simple_smoothing

bool oomph::HypreInterface::AMG_using_simple_smoothing
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 338 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(), hypre_solver_setup(), and HypreInterface().

◆ AMGEuclidSmoother_drop_tol

double oomph::HypreInterface::AMGEuclidSmoother_drop_tol

Definition at line 283 of file hypre_solver.h.

Referenced by hypre_solver_setup(), and HypreInterface().

◆ AMGEuclidSmoother_level

unsigned oomph::HypreInterface::AMGEuclidSmoother_level

Definition at line 282 of file hypre_solver.h.

Referenced by hypre_solver_setup(), and HypreInterface().

◆ AMGEuclidSmoother_print_level

unsigned oomph::HypreInterface::AMGEuclidSmoother_print_level

Definition at line 284 of file hypre_solver.h.

Referenced by hypre_solver_setup(), and HypreInterface().

◆ AMGEuclidSmoother_use_block_jacobi

bool oomph::HypreInterface::AMGEuclidSmoother_use_block_jacobi

Definition at line 279 of file hypre_solver.h.

Referenced by hypre_solver_setup(), and HypreInterface().

◆ AMGEuclidSmoother_use_ilut

bool oomph::HypreInterface::AMGEuclidSmoother_use_ilut

Definition at line 281 of file hypre_solver.h.

Referenced by hypre_solver_setup(), and HypreInterface().

◆ AMGEuclidSmoother_use_row_scaling

bool oomph::HypreInterface::AMGEuclidSmoother_use_row_scaling

Definition at line 280 of file hypre_solver.h.

Referenced by hypre_solver_setup(), and HypreInterface().

◆ Delete_input_data

bool oomph::HypreInterface::Delete_input_data
protected

Internal flag which is true when hypre_setup or hypre_solve can delete input matrix.

Definition at line 445 of file hypre_solver.h.

Referenced by hypre_matrix_setup(), oomph::HyprePreconditioner::setup(), and oomph::HypreSolver::solve().

◆ Delete_matrix

bool oomph::HypreInterface::Delete_matrix
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 465 of file hypre_solver.h.

◆ Euclid_droptol

double oomph::HypreInterface::Euclid_droptol
protected

Euclid drop tolerance for ILU(k) and ILUT factorization.

Definition at line 406 of file hypre_solver.h.

Referenced by oomph::HypreSolver::euclid_droptol(), oomph::HyprePreconditioner::euclid_droptol(), hypre_solver_setup(), and HypreInterface().

◆ Euclid_level

int oomph::HypreInterface::Euclid_level
protected

Euclid level parameter for ILU(k) factorization.

Definition at line 418 of file hypre_solver.h.

Referenced by oomph::HypreSolver::euclid_level(), oomph::HyprePreconditioner::euclid_level(), hypre_solver_setup(), and HypreInterface().

◆ Euclid_print_level

unsigned oomph::HypreInterface::Euclid_print_level
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 425 of file hypre_solver.h.

Referenced by oomph::HypreSolver::euclid_print_level(), oomph::HyprePreconditioner::euclid_print_level(), hypre_solver_setup(), and HypreInterface().

◆ Euclid_rowScale

bool oomph::HypreInterface::Euclid_rowScale
protected

◆ Euclid_using_BJ

bool oomph::HypreInterface::Euclid_using_BJ
protected

◆ Euclid_using_ILUT

bool oomph::HypreInterface::Euclid_using_ILUT
protected

◆ Existing_preconditioner

unsigned oomph::HypreInterface::Existing_preconditioner
private

Used to keep track of which preconditioner (if any) is currently stored.

Definition at line 487 of file hypre_solver.h.

Referenced by existing_preconditioner(), hypre_clean_up_memory(), hypre_solver_setup(), and HypreInterface().

◆ Existing_solver

unsigned oomph::HypreInterface::Existing_solver
private

Used to keep track of which solver (if any) is currently stored.

Definition at line 484 of file hypre_solver.h.

Referenced by existing_solver(), hypre_clean_up_memory(), hypre_solve(), hypre_solver_setup(), and HypreInterface().

◆ Hypre_distribution_pt

LinearAlgebraDistribution* oomph::HypreInterface::Hypre_distribution_pt
private

the distribution for this helpers-

Definition at line 490 of file hypre_solver.h.

Referenced by hypre_matrix_setup(), hypre_solve(), hypre_solver_setup(), HypreInterface(), and ~HypreInterface().

◆ Hypre_error_messages

bool oomph::HypreInterface::Hypre_error_messages
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:

  1. setting up the Hypre matrix
  2. setting up the preconditioner
  3. setting up the solver
  4. setting up the Hypre vectors
  5. solving
  6. getting the solution vector
  7. deallocation of solver data

Definition at line 441 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().

◆ Hypre_method

unsigned oomph::HypreInterface::Hypre_method
protected

◆ Internal_preconditioner

unsigned oomph::HypreInterface::Internal_preconditioner
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 320 of file hypre_solver.h.

Referenced by hypre_solver_setup(), HypreInterface(), oomph::HyprePreconditioner::HyprePreconditioner(), oomph::HypreSolver::internal_preconditioner(), and oomph::HyprePreconditioner::internal_preconditioner().

◆ Krylov_print_level

unsigned oomph::HypreInterface::Krylov_print_level
protected

Used to set the Hypre printing level for the Krylov subspace solvers.

Definition at line 429 of file hypre_solver.h.

Referenced by hypre_solver_setup(), HypreInterface(), and oomph::HypreSolver::krylov_print_level().

◆ Matrix_ij

HYPRE_IJMatrix oomph::HypreInterface::Matrix_ij
private

The Hypre_IJMatrix version of the matrix used in solve(...), resolve(...) or preconditioner_solve(...).

Definition at line 469 of file hypre_solver.h.

Referenced by hypre_clean_up_memory(), and hypre_matrix_setup().

◆ Matrix_par

HYPRE_ParCSRMatrix oomph::HypreInterface::Matrix_par
private

The Hypre_ParCSRMatrix version of the matrix used in solve(...), resolve(...) or preconditioner_solve(...).

Definition at line 473 of file hypre_solver.h.

Referenced by hypre_matrix_setup(), hypre_solve(), and hypre_solver_setup().

◆ Max_iter

unsigned oomph::HypreInterface::Max_iter
protected

◆ Output_info

bool oomph::HypreInterface::Output_info
protected

◆ ParaSails_filter

double oomph::HypreInterface::ParaSails_filter
protected

◆ ParaSails_nlevel

int oomph::HypreInterface::ParaSails_nlevel
protected

◆ ParaSails_symmetry

int oomph::HypreInterface::ParaSails_symmetry
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 394 of file hypre_solver.h.

Referenced by hypre_solver_setup(), HypreInterface(), oomph::HypreSolver::parasails_symmetry(), and oomph::HyprePreconditioner::parasails_symmetry().

◆ ParaSails_thresh

double oomph::HypreInterface::ParaSails_thresh
protected

◆ Preconditioner

HYPRE_Solver oomph::HypreInterface::Preconditioner
private

The internal Hypre preconditioner used in conjunction with Solver. [This is a C structure!].

Definition at line 481 of file hypre_solver.h.

◆ Returning_distributed_solution

bool oomph::HypreInterface::Returning_distributed_solution
protected

Internal flag which tell the solver if the solution Vector to be returned is distributed or not.

Definition at line 454 of file hypre_solver.h.

◆ Solver

HYPRE_Solver oomph::HypreInterface::Solver
private

The Hypre solver used in solve(...), resolve(...) or preconditioner_solve(...). [This is a C structure!].

Definition at line 477 of file hypre_solver.h.

Referenced by hypre_clean_up_memory(), hypre_solve(), and hypre_solver_setup().

◆ Tolerance

double oomph::HypreInterface::Tolerance
protected

Tolerance used to terminate solver.

Definition at line 311 of file hypre_solver.h.

Referenced by hypre_solver_setup(), HypreInterface(), oomph::HyprePreconditioner::HyprePreconditioner(), and oomph::HypreSolver::tolerance().

◆ Using_distributed_rhs

bool oomph::HypreInterface::Using_distributed_rhs
protected

Internal flag which tell the solver if the rhs Vector is distributed or not.

Definition at line 450 of file hypre_solver.h.


The documentation for this class was generated from the following files: