26 #ifndef OOMPH_HYPRE_SOLVER_HEADER
27 #define OOMPH_HYPRE_SOLVER_HEADER
30 #include "_hypre_utilities.h"
32 #include "HYPRE_parcsr_ls.h"
33 #include "HYPRE_krylov.h"
34 #include "HYPRE_IJ_mv.h"
35 #include "HYPRE_parcsr_mv.h"
51 namespace HypreHelpers
85 const LinearAlgebraDistribution* dist_pt,
86 HYPRE_IJVector& hypre_ij_vector,
87 HYPRE_ParVector& hypre_par_vector);
97 HYPRE_IJVector& hypre_ij_vector,
98 HYPRE_ParVector& hypre_par_vector);
103 HYPRE_IJMatrix& hypre_ij_matrix,
104 HYPRE_ParCSRMatrix& hypre_par_matrix,
105 LinearAlgebraDistribution* dist_pt);
110 const bool& use_row_scaling,
111 const bool& use_ilut,
113 const double& drop_tol,
114 const int& print_level,
115 HYPRE_Solver& euclid_object);
140 #ifndef HYPRE_SEQUENTIAL
144 std::ostringstream error_message;
145 error_message <<
"When using the MPI version of Hypre please first\n"
146 <<
"call function MPI_Helpers::setup()\n";
148 OOMPH_CURRENT_FUNCTION,
149 OOMPH_EXCEPTION_LOCATION);
867 <<
"\": My_cumulative_preconditioner_solve_time = "
903 static std::map<std::string, unsigned>
1245 namespace Hypre_default_settings
A class for compressed row matrices. This is a distributable object.
Abstract base class for matrices of doubles – adds abstract interfaces for solving,...
A vector in the mathematical sense, initially developed for linear algebra type applications....
An interface class to the suite of Hypre solvers and preconditioners to allow use of:
bool Output_info
Flag is true to output info and results of timings.
bool Delete_input_data
Internal flag which is true when hypre_setup or hypre_solve can delete input matrix.
unsigned AMG_coarsening
AMG coarsening strategy. Coarsening types include: 0 = CLJP (parallel coarsening using independent se...
int Euclid_level
Euclid level parameter for ILU(k) factorization.
double AMGEuclidSmoother_drop_tol
bool Euclid_using_BJ
Flag to determine if Block Jacobi is used instead of PILU.
bool Euclid_rowScale
Flag to switch on Euclid row scaling.
bool Returning_distributed_solution
Internal flag which tell the solver if the solution Vector to be returned is distributed or not.
double AMG_strength
Connection strength threshold parameter for BoomerAMG.
HypreInterface(const HypreInterface &)=delete
Broken copy constructor.
void enable_hypre_error_messages()
Turn on the hypre error messages.
LinearAlgebraDistribution * Hypre_distribution_pt
the distribution for this helpers-
int ParaSails_nlevel
ParaSails nlevel parameter.
double ParaSails_filter
ParaSails filter parameter.
bool AMGEuclidSmoother_use_row_scaling
unsigned AMGEuclidSmoother_level
HYPRE_Solver Solver
The Hypre solver used in solve(...), resolve(...) or preconditioner_solve(...). [This is a C structur...
int ParaSails_symmetry
ParaSails symmetry flag, used to inform ParaSails of Symmetry of definitenss of problem and type of P...
bool Using_distributed_rhs
Internal flag which tell the solver if the rhs Vector is distributed or not.
HypreInterface()
Constructor.
unsigned existing_preconditioner()
Function return value of which preconditioner (if any) is stored.
HYPRE_ParCSRMatrix Matrix_par
The Hypre_ParCSRMatrix version of the matrix used in solve(...), resolve(...) or preconditioner_solve...
bool AMGEuclidSmoother_use_ilut
void operator=(const HypreInterface &)=delete
Broken assignment operator.
void disable_hypre_error_messages()
Turn off hypre error messages.
unsigned Internal_preconditioner
Preconditioner method flag used with Hypre's PCG, GMRES and BiCGStab in solve(...) or resolve(....
void hypre_solve(const DoubleVector &rhs, DoubleVector &solution)
Helper function performs a solve if any solver exists.
Hypre_method_types
Enumerated flag to define which Hypre methods are used CAREFUL: DON'T CHANGE THE ORDER OF THESE!
~HypreInterface()
Destructor.
HYPRE_Solver Preconditioner
The internal Hypre preconditioner used in conjunction with Solver. [This is a C structure!...
unsigned AMGEuclidSmoother_print_level
bool AMGEuclidSmoother_use_block_jacobi
unsigned AMG_complex_smoother
Complex smoothing methods used in BoomerAMG. Relaxation types are: 6 = Schwarz 7 = Pilut 8 = ParaSail...
bool Delete_matrix
Hypre copies matrix data from oomph-lib's CRDoubleMatrix into its own data structures,...
double ParaSails_thresh
ParaSails thresh parameter.
double AMG_truncation
Interpolation truncation factor for BoomerAMG.
unsigned Krylov_print_level
Used to set the Hypre printing level for the Krylov subspace solvers.
unsigned AMG_smoother_iterations
The number of smoother iterations to apply.
unsigned existing_solver()
Function to return value of which solver (if any) is currently stored.
void hypre_clean_up_memory()
Function deletes all solver data.
bool Hypre_error_messages
Flag to determine if non-zero values of the Hypre error flag plus Hypre error messages are output to ...
double Tolerance
Tolerance used to terminate solver.
bool AMG_using_simple_smoothing
Flag to determine whether simple smoothers (determined by the AMG_simple_smoother flag) or complex sm...
unsigned Existing_solver
Used to keep track of which solver (if any) is currently stored.
unsigned Max_iter
Maximum number of iterations used in solver.
double AMG_damping
Damping factor for BoomerAMG smoothed Jacobi or hybrid SOR.
unsigned AMG_simple_smoother
Simple smoothing methods used in BoomerAMG. Relaxation types include: 0 = Jacobi 1 = Gauss-Seidel,...
unsigned AMG_print_level
Used to set the Hypre printing level for AMG 0: no printout 1: print setup information 2: print solve...
HYPRE_IJMatrix Matrix_ij
The Hypre_IJMatrix version of the matrix used in solve(...), resolve(...) or preconditioner_solve(....
unsigned Euclid_print_level
Flag to set the level of printing from Euclid when the Euclid destructor is called 0: no printing (de...
unsigned AMG_max_levels
Maximum number of levels used in AMG.
unsigned Hypre_method
Hypre method flag. Valid values are specified in enumeration.
void hypre_matrix_setup(CRDoubleMatrix *matrix_pt)
Function which sets values of First_global_row, Last_global_row and other partitioning data and creat...
double AMG_max_row_sum
Parameter to identify diagonally dominant parts of the matrix in AMG.
double Euclid_droptol
Euclid drop tolerance for ILU(k) and ILUT factorization.
unsigned Existing_preconditioner
Used to keep track of which preconditioner (if any) is currently stored.
bool Euclid_using_ILUT
Flag to determine if ILUT (if true) or ILU(k) is used in Euclid.
void hypre_solver_setup()
Sets up the data required for to use as an oomph-lib LinearSolver or Preconditioner....
///////////////////////////////////////////////////////////////////// ///////////////////////////////...
~HyprePreconditioner()
Destructor.
unsigned & amg_smoother_iterations()
Access function to AMG_smoother_iterations.
unsigned & amg_simple_smoother()
Access function to AMG_simple_smoother flag.
double & amg_max_row_sum()
Access function to AMG_max_row_sum.
static void reset_cumulative_solve_times()
Reset cumulative solve times.
double & amg_strength()
Access function to AMG_strength.
static std::map< std::string, unsigned > Context_based_cumulative_npreconditioner_solve
Static unsigned that accumulates the number of preconditioner solves of all instantiations of this cl...
unsigned & amg_iterations()
Return function for Max_iter.
void disable_euclid_using_BJ()
Disable use of Block Jacobi, so PILU will be used.
std::string Context_string
String can be used to provide context for annotation.
double & parasails_thresh()
Access function to ParaSails thresh parameter.
int & parasails_nlevel()
Access function to ParaSails nlevel parameter.
void enable_delete_matrix()
Hypre copies matrix data from oomph-lib's CRDoubleMatrix or DistributedCRDoubleMatrix into its own da...
void disable_report_my_cumulative_preconditioner_solve_time()
Disable reporting of cumulative solve time in destructor.
double & amg_damping()
Access function to AMG_damping parameter.
unsigned & amg_print_level()
Access function to AMG_print_level.
void set_amg_iterations(const unsigned &amg_iterations)
Function to set the number of times to apply BoomerAMG.
void clean_up_memory()
Function deletes all solver data.
void disable_delete_matrix()
Hypre copies matrix data from oomph-lib's CRDoubleMatrix or DistributedCRDoubleMatrix into its own da...
static double Cumulative_preconditioner_solve_time
Static double that accumulates the preconditioner solve time of all instantiations of this class....
void euclid_using_ILUT()
Function to switch on ILUT factorization for Euclid (default is ILU(k) factorization)
void setup()
Function to set up a preconditioner for the linear system defined by matrix_pt. This function is requ...
void enable_euclid_rowScale()
Enable euclid rowScaling.
unsigned & hypre_method()
Access function to Hypre_method flag – specified via enumeration.
void euclid_using_ILUK()
Function to switch on ILU(k) factorization for Euclid (default is ILU(k) factorization)
int & euclid_level()
Access function to Euclid level parameter for ILU(k) factorization.
bool Report_my_cumulative_preconditioner_solve_time
Bool to request report of My_cumulative_preconditioner_solve_time in destructor.
void amg_using_complex_smoothing()
Function to select use of 'complex' AMG smoothers as controlled by the flag AMG_complex_smoother.
static std::map< std::string, double > Context_based_cumulative_solve_time
Static double that accumulates the preconditioner solve time of all instantiations of this class,...
unsigned & amg_max_levels()
Access function to AMG_max_levels.
void enable_euclid_using_BJ()
Enable use of Block Jacobi as opposed to PILU.
double & amg_truncation()
Access function to AMG_truncation.
bool Delete_matrix
Hypre copies matrix data from oomph-lib's CRDoubleMatrix or DistributedCRDoubleMatrix into its own da...
void enable_doc_time()
Enable documentation of preconditioner timings.
void disable_doc_time()
Disable documentation of preconditioner timings.
static std::map< std::string, unsigned > Context_based_nrow
Static unsigned that stores nrow for the most recent instantiations of this class,...
void amg_using_simple_smoothing()
Function to select use of 'simple' AMG smoothers as controlled by the flag AMG_simple_smoother.
bool & amg_using_simple_smoothing_flag()
Return function for the AMG_using_simple_smoothing_flag.
double My_cumulative_preconditioner_solve_time
Private double that accumulates the preconditioner solve time of thi instantiation of this class....
void use_Euclid()
Function to select use Euclid as the preconditioner.
double & euclid_droptol()
Access function to Euclid drop tolerance parameter.
unsigned & euclid_print_level()
Function to set the level of printing from Euclid when the Euclid destructor is called 0: no printing...
void disable_euclid_rowScale()
Disable euclid row scaling.
void operator=(const HyprePreconditioner &)=delete
Broken assignment operator.
void use_BoomerAMG()
Function to select BoomerAMG as the preconditioner.
unsigned & internal_preconditioner()
Access function to Internal_preconditioner flag – specified via enumeration.
static void report_cumulative_solve_times()
Report cumulative solve times of all instantiations of this class.
void use_ParaSails()
Function to select ParaSails as the preconditioner.
int & parasails_symmetry()
Access function to ParaSails symmetry flag.
unsigned & amg_complex_smoother()
Access function to AMG_complex_smoother flag.
unsigned & amg_coarsening()
Access function to AMG_coarsening flag.
static unsigned Cumulative_npreconditioner_solve
Static unsigned that accumulates the number of preconditioner solves of all instantiations of this cl...
void preconditioner_solve(const DoubleVector &r, DoubleVector &z)
Function applies solver to vector r for preconditioning. This requires a call to setup(....
HyprePreconditioner(const HyprePreconditioner &)=delete
Broken copy constructor.
HyprePreconditioner(const std::string &context_string="")
Constructor. Provide optional string that is used in annotation of performance.
double & parasails_filter()
Access function to ParaSails filter parameter.
void enable_report_my_cumulative_preconditioner_solve_time()
Enable reporting of cumulative solve time in destructor.
///////////////////////////////////////////////////////////////////// ///////////////////////////////...
unsigned & amg_simple_smoother()
Access function to AMG_simple_smoother flag.
void disable_euclid_using_BJ()
Disable use of Block Jacobi, so PILU will be used.
int & euclid_level()
Access function to Euclid level parameter for ILU(k) factorization.
void disable_euclid_rowScale()
Disable euclid row scaling.
double & amg_damping()
Access function to AMG_damping parameter.
unsigned & max_iter()
Access function to Max_iter.
void enable_euclid_using_BJ()
Enable use of Block Jacobi as opposed to PILU.
double & euclid_droptol()
Access function to Euclid drop tolerance parameter.
unsigned & hypre_method()
Access function to Hypre_method flag – specified via enumeration.
unsigned & euclid_print_level()
Function to set the level of printing from Euclid when the Euclid destructor is called 0: no printing...
void solve(Problem *const &problem_pt, DoubleVector &solution)
Function which uses problem_pt's get_jacobian(...) function to generate a linear system which is then...
double & amg_max_row_sum()
Access function to AMG_max_row_sum.
unsigned & amg_print_level()
Access function to AMG_print_level.
HypreSolver(const HypreSolver &)=delete
Broken copy constructor.
~HypreSolver()
Empty destructor.
unsigned & krylov_print_level()
Access function to Krylov_print_level.
double & tolerance()
Access function to Tolerance.
void clean_up_memory()
Function deletes all solver data.
unsigned & amg_coarsening()
Access function to AMG_coarsening flag.
double & parasails_thresh()
Access function to ParaSails thresh parameter.
HypreSolver()
Constructor.
double & parasails_filter()
Access function to ParaSails filter parameter.
unsigned & amg_complex_smoother()
Access function to AMG_complex_smoother flag.
bool Delete_matrix
Hypre copies matrix data from oomph-lib's CRDoubleMatrix or DistributedCRDoubleMatrix into its own da...
void amg_using_complex_smoothing()
Function to select use of 'complex' AMG smoothers as controlled by AMG_complex_smoother flag.
void euclid_using_ILUT()
Function to switch on ILUT factorization for Euclid (default is ILU(k) factorization)
void amg_using_simple_smoothing()
Function to select use of 'simple' AMG smoothers as controlled by AMG_simple_smoother flag.
void resolve(const DoubleVector &rhs, DoubleVector &solution)
Function to resolve a linear system using the existing solver data, allowing a solve with a new right...
void enable_euclid_rowScale()
Enable euclid rowScaling.
unsigned & amg_max_levels()
Access function to AMG_max_levels.
int & parasails_symmetry()
Access function to ParaSails symmetry flag.
unsigned & amg_smoother_iterations()
Access function to AMG_smoother_iterations.
unsigned & internal_preconditioner()
Access function to Internal_preconditioner flag – specified via enumeration.
void euclid_using_ILUK()
Function to switch on ILU(k) factorization for Euclid (default is ILU(k) factorization)
void operator=(const HypreSolver &)=delete
Broken assignment operator.
double & amg_truncation()
Access function to AMG_truncation.
int & parasails_nlevel()
Access function to ParaSails nlevel parameter.
void disable_resolve()
Disable resolve function (overloads the LinearSolver disable_resolve function).
void enable_delete_matrix()
Hypre copies matrix data from oomph-lib's CRDoubleMatrix or DistributedCRDoubleMatrix into its own da...
void disable_delete_matrix()
Hypre copies matrix data from oomph-lib's CRDoubleMatrix or DistributedCRDoubleMatrix into its own da...
double & amg_strength()
Access function to AMG_strength.
Describes the distribution of a distributable linear algebra type object. Typically this is a contain...
Base class for all linear solvers. This merely defines standard interfaces for linear solvers,...
bool Doc_time
Boolean flag that indicates whether the time taken.
bool Enable_resolve
Boolean that indicates whether the matrix (or its factors, in the case of direct solver) should be st...
static bool mpi_has_been_initialised()
return true if MPI has been initialised
static OomphCommunicator * communicator_pt()
access to global communicator. This is the oomph-lib equivalent of MPI_COMM_WORLD
An OomphLibError object which should be thrown when an run-time error is encountered....
Preconditioner base class. Gives an interface to call all other preconditioners through and stores th...
////////////////////////////////////////////////////////////////// //////////////////////////////////...
std::string string(const unsigned &i)
Return the i-th string or "" if the relevant string hasn't been defined.
void create_HYPRE_Matrix(CRDoubleMatrix *oomph_matrix, HYPRE_IJMatrix &hypre_ij_matrix, HYPRE_ParCSRMatrix &hypre_par_matrix, LinearAlgebraDistribution *dist_pt)
Helper function to create a serial HYPRE_IJMatrix and HYPRE_ParCSRMatrix from a CRDoubleMatrix NOTE: ...
void create_HYPRE_Vector(const DoubleVector &oomph_vec, const LinearAlgebraDistribution *dist_pt, HYPRE_IJVector &hypre_ij_vector, HYPRE_ParVector &hypre_par_vector)
Helper function to create a HYPRE_IJVector and HYPRE_ParVector.
int check_HYPRE_error_flag(std::ostringstream &message)
Helper function to check the Hypre error flag, return the message associated with any error,...
unsigned AMG_coarsening
Default AMG coarsening strategy. Coarsening types include: 0 = CLJP (parallel coarsening using indepe...
double AMG_truncation
AMG interpolation truncation factor.
double AMG_strength
Default for AMG strength (0.25 recommended for 2D problems; larger (0.5-0.75, say) for 3D.
void euclid_settings_helper(const bool &use_block_jacobi, const bool &use_row_scaling, const bool &use_ilut, const int &level, const double &drop_tol, const int &print_level, HYPRE_Solver &euclid_object)
Helper function to set Euclid options using a command line like array.
void set_defaults_for_navier_stokes_momentum_block(HyprePreconditioner *hypre_preconditioner_pt)
Set default parameters for use as preconditioner in for momentum block in Navier-Stokes problem.
void set_defaults_for_3D_poisson_problem(HyprePreconditioner *hypre_preconditioner_pt)
Set default parameters for use as preconditioner in 3D Poisson-type problem.
void set_defaults_for_2D_poisson_problem(HyprePreconditioner *hypre_preconditioner_pt)
Set default parameters for use as preconditioner in 2D Poisson-type problem.
//////////////////////////////////////////////////////////////////// ////////////////////////////////...
OomphInfo oomph_info
Single (global) instantiation of the OomphInfo object – this is used throughout the library as a "rep...