///////////////////////////////////////////////////////////////////// ///////////////////////////////////////////////////////////////////// More...
#include <hypre_solver.h>
 Inheritance diagram for oomph::HyprePreconditioner:
 Inheritance diagram for oomph::HyprePreconditioner:| Public Member Functions | |
| HyprePreconditioner (const std::string &context_string="") | |
| Constructor. Provide optional string that is used in annotation of performance.  More... | |
| ~HyprePreconditioner () | |
| Destructor.  More... | |
| HyprePreconditioner (const HyprePreconditioner &)=delete | |
| Broken copy constructor.  More... | |
| void | operator= (const HyprePreconditioner &)=delete | 
| Broken assignment operator.  More... | |
| void | enable_report_my_cumulative_preconditioner_solve_time () | 
| Enable reporting of cumulative solve time in destructor.  More... | |
| void | disable_report_my_cumulative_preconditioner_solve_time () | 
| Disable reporting of cumulative solve time in destructor.  More... | |
| void | enable_doc_time () | 
| Enable documentation of preconditioner timings.  More... | |
| void | disable_doc_time () | 
| Disable documentation of preconditioner timings.  More... | |
| unsigned & | hypre_method () | 
| Access function to Hypre_method flag – specified via enumeration.  More... | |
| unsigned & | internal_preconditioner () | 
| Access function to Internal_preconditioner flag – specified via enumeration.  More... | |
| void | use_BoomerAMG () | 
| Function to select BoomerAMG as the preconditioner.  More... | |
| void | set_amg_iterations (const unsigned &amg_iterations) | 
| Function to set the number of times to apply BoomerAMG.  More... | |
| unsigned & | amg_iterations () | 
| Return function for Max_iter.  More... | |
| void | amg_using_simple_smoothing () | 
| Function to select use of 'simple' AMG smoothers as controlled by the flag AMG_simple_smoother.  More... | |
| unsigned & | amg_simple_smoother () | 
| Access function to AMG_simple_smoother flag.  More... | |
| void | amg_using_complex_smoothing () | 
| Function to select use of 'complex' AMG smoothers as controlled by the flag AMG_complex_smoother.  More... | |
| unsigned & | amg_complex_smoother () | 
| Access function to AMG_complex_smoother flag.  More... | |
| bool & | amg_using_simple_smoothing_flag () | 
| Return function for the AMG_using_simple_smoothing_flag.  More... | |
| unsigned & | amg_print_level () | 
| Access function to AMG_print_level.  More... | |
| unsigned & | amg_smoother_iterations () | 
| Access function to AMG_smoother_iterations.  More... | |
| unsigned & | amg_coarsening () | 
| Access function to AMG_coarsening flag.  More... | |
| unsigned & | amg_max_levels () | 
| Access function to AMG_max_levels.  More... | |
| double & | amg_damping () | 
| Access function to AMG_damping parameter.  More... | |
| double & | amg_strength () | 
| Access function to AMG_strength.  More... | |
| double & | amg_max_row_sum () | 
| Access function to AMG_max_row_sum.  More... | |
| double & | amg_truncation () | 
| Access function to AMG_truncation.  More... | |
| void | use_ParaSails () | 
| Function to select ParaSails as the preconditioner.  More... | |
| int & | parasails_symmetry () | 
| Access function to ParaSails symmetry flag.  More... | |
| int & | parasails_nlevel () | 
| Access function to ParaSails nlevel parameter.  More... | |
| double & | parasails_thresh () | 
| Access function to ParaSails thresh parameter.  More... | |
| double & | parasails_filter () | 
| Access function to ParaSails filter parameter.  More... | |
| void | use_Euclid () | 
| Function to select use Euclid as the preconditioner.  More... | |
| double & | euclid_droptol () | 
| Access function to Euclid drop tolerance parameter.  More... | |
| int & | euclid_level () | 
| Access function to Euclid level parameter for ILU(k) factorization.  More... | |
| void | enable_euclid_rowScale () | 
| Enable euclid rowScaling.  More... | |
| void | disable_euclid_rowScale () | 
| Disable euclid row scaling.  More... | |
| void | enable_euclid_using_BJ () | 
| Enable use of Block Jacobi as opposed to PILU.  More... | |
| void | disable_euclid_using_BJ () | 
| Disable use of Block Jacobi, so PILU will be used.  More... | |
| void | euclid_using_ILUK () | 
| Function to switch on ILU(k) factorization for Euclid (default is ILU(k) factorization)  More... | |
| void | euclid_using_ILUT () | 
| Function to switch on ILUT factorization for Euclid (default is ILU(k) factorization)  More... | |
| unsigned & | euclid_print_level () | 
| Function 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... | |
| void | enable_delete_matrix () | 
| Hypre copies matrix data from oomph-lib's CRDoubleMatrix or DistributedCRDoubleMatrix into its own data structures, doubling the memory requirements for the matrix. 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 calling this function which changes the flag from false (its default) to true.  More... | |
| void | disable_delete_matrix () | 
| Hypre copies matrix data from oomph-lib's CRDoubleMatrix or DistributedCRDoubleMatrix into its own data structures, doubling the memory requirements for the matrix. Calling this function ensures that the matrix is not deleted.  More... | |
| void | setup () | 
| Function to set up a preconditioner for the linear system defined by matrix_pt. This function is required when preconditioning and must be called before using the preconditioner_solve(...) function. This interface allows HyprePreconditioner to be used as a Preconditioner object for oomph-lib's own IterativeLinearSolver class. Note: matrix_pt must point to an object of type CRDoubleMatrix or DistributedCRDoubleMatrix. Note: Hypre copies matrix data from oomph-lib's CRDoubleMatrix and DistributedCRDoubleMatrix into its own data structures, doubling the memory requirements for the matrix. 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 calling the enable_delete_matrix() function.  More... | |
| void | preconditioner_solve (const DoubleVector &r, DoubleVector &z) | 
| Function applies solver to vector r for preconditioning. This requires a call to setup(...) first. Note: Hypre copies matrix data from oomph-lib's CRDoubleMatrix or DistributedCRDoubleMatrix into its own data structures, doubling the memory requirements for the matrix. 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 with the delete_matrix() function.  More... | |
| void | clean_up_memory () | 
| Function deletes all solver data.  More... | |
|  Public Member Functions inherited from oomph::Preconditioner | |
| Preconditioner () | |
| Constructor.  More... | |
| Preconditioner (const Preconditioner &)=delete | |
| Broken copy constructor.  More... | |
| void | operator= (const Preconditioner &)=delete | 
| Broken assignment operator.  More... | |
| virtual | ~Preconditioner () | 
| Destructor (empty)  More... | |
| virtual void | preconditioner_solve_transpose (const DoubleVector &r, DoubleVector &z) | 
| Apply the preconditioner. Pure virtual generic interface function. This method should apply the preconditioner operator to the vector r and return the vector z. (broken virtual)  More... | |
| void | setup (DoubleMatrixBase *matrix_pt) | 
| Setup the preconditioner: store the matrix pointer and the communicator pointer then call preconditioner specific setup() function.  More... | |
| void | setup (const Problem *problem_pt, DoubleMatrixBase *matrix_pt) | 
| Compatability layer for old preconditioners where problem pointers were needed. The problem pointer is only used to get a communicator pointer.  More... | |
| void | enable_silent_preconditioner_setup () | 
| Set up the block preconditioner quietly!  More... | |
| void | disable_silent_preconditioner_setup () | 
| Be verbose in the block preconditioner setup.  More... | |
| virtual DoubleMatrixBase * | matrix_pt () const | 
| Get function for matrix pointer.  More... | |
| virtual void | set_matrix_pt (DoubleMatrixBase *matrix_pt) | 
| Set the matrix pointer.  More... | |
| virtual const OomphCommunicator * | comm_pt () const | 
| Get function for comm pointer.  More... | |
| virtual void | set_comm_pt (const OomphCommunicator *const comm_pt) | 
| Set the communicator pointer.  More... | |
| double | setup_time () const | 
| Returns the time to setup the preconditioner.  More... | |
| virtual void | turn_into_subsidiary_block_preconditioner (BlockPreconditioner< CRDoubleMatrix > *master_block_prec_pt, const Vector< unsigned > &doftype_in_master_preconditioner_coarse) | 
| Virtual interface function for making a preconditioner a subsidiary of a block preconditioner. By default nothing is needed, but if this preconditioner is also a block preconditioner then things need to happen. There's an assumption here that the block preconditioner will be in CR form but since that assumption is hard coded all over BlockPreconditioner we're safe.  More... | |
| virtual void | turn_into_subsidiary_block_preconditioner (BlockPreconditioner< CRDoubleMatrix > *master_block_prec_pt, const Vector< unsigned > &doftype_in_master_preconditioner_coarse, const Vector< Vector< unsigned >> &doftype_coarsen_map_coarse) | 
| Virtual interface function for making a preconditioner a subsidiary of a block preconditioner. By default nothing is needed, but if this preconditioner is also a block preconditioner then things need to happen. Version for coarsening dof-types.  More... | |
|  Public Member Functions inherited from oomph::DistributableLinearAlgebraObject | |
| DistributableLinearAlgebraObject () | |
| Default constructor - create a distribution.  More... | |
| DistributableLinearAlgebraObject (const DistributableLinearAlgebraObject &matrix)=delete | |
| Broken copy constructor.  More... | |
| void | operator= (const DistributableLinearAlgebraObject &)=delete | 
| Broken assignment operator.  More... | |
| virtual | ~DistributableLinearAlgebraObject () | 
| Destructor.  More... | |
| LinearAlgebraDistribution * | distribution_pt () const | 
| access to the LinearAlgebraDistribution  More... | |
| unsigned | nrow () const | 
| access function to the number of global rows.  More... | |
| unsigned | nrow_local () const | 
| access function for the num of local rows on this processor.  More... | |
| unsigned | nrow_local (const unsigned &p) const | 
| access function for the num of local rows on this processor.  More... | |
| unsigned | first_row () const | 
| access function for the first row on this processor  More... | |
| unsigned | first_row (const unsigned &p) const | 
| access function for the first row on this processor  More... | |
| bool | distributed () const | 
| distribution is serial or distributed  More... | |
| bool | distribution_built () const | 
| if the communicator_pt is null then the distribution is not setup then false is returned, otherwise return true  More... | |
| void | build_distribution (const LinearAlgebraDistribution *const dist_pt) | 
| setup the distribution of this distributable linear algebra object  More... | |
| void | build_distribution (const LinearAlgebraDistribution &dist) | 
| setup the distribution of this distributable linear algebra object  More... | |
|  Public Member Functions inherited from oomph::HypreInterface | |
| 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... | |
| Static Public Member Functions | |
| static void | report_cumulative_solve_times () | 
| Report cumulative solve times of all instantiations of this class.  More... | |
| static void | reset_cumulative_solve_times () | 
| Reset cumulative solve times.  More... | |
| Static Public Attributes | |
| static double | Cumulative_preconditioner_solve_time = 0.0 | 
| Static double that accumulates the preconditioner solve time of all instantiations of this class. Reset it manually, e.g. after every Newton solve, using reset_cumulative_solve_times().  More... | |
| 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, labeled by context string. Reset it manually, e.g. after every Newton solve, using reset_cumulative_solve_times().  More... | |
| static unsigned | Cumulative_npreconditioner_solve = 0 | 
| Static unsigned that accumulates the number of preconditioner solves of all instantiations of this class. Reset it manually, e.g. after every Newton solve, using reset_cumulative_solve_times().  More... | |
| 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 class, labeled by context string. Reset it manually, e.g. after every Newton solve, using reset_cumulative_solve_times().  More... | |
| static std::map< std::string, unsigned > | Context_based_nrow | 
| Static unsigned that stores nrow for the most recent instantiations of this class, labeled by context string. Reset it manually, e.g. after every Newton solve, using reset_cumulative_solve_times().  More... | |
| Private Attributes | |
| bool | Delete_matrix | 
| Hypre copies matrix data from oomph-lib's CRDoubleMatrix or DistributedCRDoubleMatrix into its own data structures, doubling the memory requirements for the matrix. 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... | |
| bool | Doc_time | 
| double | My_cumulative_preconditioner_solve_time | 
| Private double that accumulates the preconditioner solve time of thi instantiation of this class. Is reported in destructor if Report_my_cumulative_preconditioner_solve_time is set to true.  More... | |
| bool | Report_my_cumulative_preconditioner_solve_time | 
| Bool to request report of My_cumulative_preconditioner_solve_time in destructor.  More... | |
| std::string | Context_string | 
| String can be used to provide context for annotation.  More... | |
| Additional Inherited Members | |
|  Public Types inherited from oomph::HypreInterface | |
| 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 Attributes inherited from oomph::HypreInterface | |
| 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 inherited from oomph::DistributableLinearAlgebraObject | |
| void | clear_distribution () | 
| clear the distribution of this distributable linear algebra object  More... | |
|  Protected Member Functions inherited from oomph::HypreInterface | |
| 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 inherited from oomph::Preconditioner | |
| bool | Silent_preconditioner_setup | 
| Boolean to indicate whether or not the build should be done silently.  More... | |
| std::ostream * | Stream_pt | 
| Pointer to the output stream – defaults to std::cout.  More... | |
|  Protected Attributes inherited from oomph::HypreInterface | |
| 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... | |
///////////////////////////////////////////////////////////////////// /////////////////////////////////////////////////////////////////////
An Preconditioner class using the suite of Hypre preconditioners to allow
BoomerAMG (Algebraic Multi Grid), Euclid (ILU) or ParaSails (Approximate inverse)
to be used for Preconditioner::preconditioner_solve(...). By default BoomerAMG is used.
Definition at line 825 of file hypre_solver.h.
| 
 | inline | 
Constructor. Provide optional string that is used in annotation of performance.
Definition at line 830 of file hypre_solver.h.
References oomph::HypreInterface::BoomerAMG, Context_string, Delete_matrix, Doc_time, oomph::HypreInterface::Hypre_method, oomph::HypreInterface::Internal_preconditioner, oomph::HypreInterface::Max_iter, My_cumulative_preconditioner_solve_time, oomph::HypreInterface::None, Report_my_cumulative_preconditioner_solve_time, and oomph::HypreInterface::Tolerance.
| 
 | inline | 
Destructor.
Definition at line 862 of file hypre_solver.h.
References Context_string, My_cumulative_preconditioner_solve_time, oomph::oomph_info, and Report_my_cumulative_preconditioner_solve_time.
| 
 | delete | 
Broken copy constructor.
| 
 | inline | 
Access function to AMG_coarsening flag.
Definition at line 1020 of file hypre_solver.h.
References oomph::HypreInterface::AMG_coarsening.
Referenced by oomph::Lagrange_Enforced_Flow_Preconditioner_Subsidiary_Operator_Helper::boomer_amg2v22_for_3D_momentum(), oomph::Lagrange_Enforced_Flow_Preconditioner_Subsidiary_Operator_Helper::boomer_amg_for_2D_momentum_simple_visc(), oomph::Lagrange_Enforced_Flow_Preconditioner_Subsidiary_Operator_Helper::boomer_amg_for_2D_momentum_stressdiv_visc(), oomph::Lagrange_Enforced_Flow_Preconditioner_Subsidiary_Operator_Helper::boomer_amg_for_2D_poisson_problem(), oomph::Lagrange_Enforced_Flow_Preconditioner_Subsidiary_Operator_Helper::boomer_amg_for_3D_momentum(), oomph::Lagrange_Enforced_Flow_Preconditioner_Subsidiary_Operator_Helper::boomer_amg_for_3D_poisson_problem(), oomph::Pseudo_Elastic_Preconditioner_Subsidiary_Operator_Helper::get_elastic_preconditioner_hypre(), oomph::Biharmonic_schur_complement_Hypre_defaults::set_defaults(), and oomph::Hypre_default_settings::set_defaults_for_2D_poisson_problem().
| 
 | inline | 
Access function to AMG_complex_smoother flag.
Definition at line 996 of file hypre_solver.h.
References oomph::HypreInterface::AMG_complex_smoother.
| 
 | inline | 
Access function to AMG_damping parameter.
Definition at line 1032 of file hypre_solver.h.
References oomph::HypreInterface::AMG_damping.
Referenced by oomph::Lagrange_Enforced_Flow_Preconditioner_Subsidiary_Operator_Helper::boomer_amg2v22_for_3D_momentum(), oomph::Lagrange_Enforced_Flow_Preconditioner_Subsidiary_Operator_Helper::boomer_amg_for_2D_momentum_simple_visc(), oomph::Lagrange_Enforced_Flow_Preconditioner_Subsidiary_Operator_Helper::boomer_amg_for_2D_momentum_stressdiv_visc(), oomph::Lagrange_Enforced_Flow_Preconditioner_Subsidiary_Operator_Helper::boomer_amg_for_2D_poisson_problem(), oomph::Lagrange_Enforced_Flow_Preconditioner_Subsidiary_Operator_Helper::boomer_amg_for_3D_momentum(), oomph::Lagrange_Enforced_Flow_Preconditioner_Subsidiary_Operator_Helper::boomer_amg_for_3D_poisson_problem(), oomph::Pseudo_Elastic_Preconditioner_Subsidiary_Operator_Helper::get_elastic_preconditioner_hypre(), and oomph::Hypre_default_settings::set_defaults_for_navier_stokes_momentum_block().
| 
 | inline | 
Return function for Max_iter.
Definition at line 970 of file hypre_solver.h.
References oomph::HypreInterface::Max_iter.
Referenced by set_amg_iterations(), and oomph::Biharmonic_schur_complement_Hypre_defaults::set_defaults().
| 
 | inline | 
Access function to AMG_max_levels.
Definition at line 1026 of file hypre_solver.h.
References oomph::HypreInterface::AMG_max_levels.
| 
 | inline | 
Access function to AMG_max_row_sum.
Definition at line 1044 of file hypre_solver.h.
References oomph::HypreInterface::AMG_max_row_sum.
| 
 | inline | 
Access function to AMG_print_level.
Definition at line 1008 of file hypre_solver.h.
References oomph::HypreInterface::AMG_print_level.
| 
 | inline | 
Access function to AMG_simple_smoother flag.
Definition at line 983 of file hypre_solver.h.
References oomph::HypreInterface::AMG_simple_smoother.
Referenced by oomph::Lagrange_Enforced_Flow_Preconditioner_Subsidiary_Operator_Helper::boomer_amg2v22_for_3D_momentum(), oomph::Lagrange_Enforced_Flow_Preconditioner_Subsidiary_Operator_Helper::boomer_amg_for_2D_momentum_simple_visc(), oomph::Lagrange_Enforced_Flow_Preconditioner_Subsidiary_Operator_Helper::boomer_amg_for_2D_momentum_stressdiv_visc(), oomph::Lagrange_Enforced_Flow_Preconditioner_Subsidiary_Operator_Helper::boomer_amg_for_2D_poisson_problem(), oomph::Lagrange_Enforced_Flow_Preconditioner_Subsidiary_Operator_Helper::boomer_amg_for_3D_momentum(), oomph::Lagrange_Enforced_Flow_Preconditioner_Subsidiary_Operator_Helper::boomer_amg_for_3D_poisson_problem(), oomph::Pseudo_Elastic_Preconditioner_Subsidiary_Operator_Helper::get_elastic_preconditioner_hypre(), oomph::Biharmonic_schur_complement_Hypre_defaults::set_defaults(), oomph::Hypre_default_settings::set_defaults_for_2D_poisson_problem(), and oomph::Hypre_default_settings::set_defaults_for_navier_stokes_momentum_block().
| 
 | inline | 
Access function to AMG_smoother_iterations.
Definition at line 1014 of file hypre_solver.h.
References oomph::HypreInterface::AMG_smoother_iterations.
Referenced by oomph::Lagrange_Enforced_Flow_Preconditioner_Subsidiary_Operator_Helper::boomer_amg2v22_for_3D_momentum(), oomph::Lagrange_Enforced_Flow_Preconditioner_Subsidiary_Operator_Helper::boomer_amg_for_2D_momentum_simple_visc(), oomph::Lagrange_Enforced_Flow_Preconditioner_Subsidiary_Operator_Helper::boomer_amg_for_2D_momentum_stressdiv_visc(), oomph::Lagrange_Enforced_Flow_Preconditioner_Subsidiary_Operator_Helper::boomer_amg_for_2D_poisson_problem(), oomph::Lagrange_Enforced_Flow_Preconditioner_Subsidiary_Operator_Helper::boomer_amg_for_3D_momentum(), oomph::Lagrange_Enforced_Flow_Preconditioner_Subsidiary_Operator_Helper::boomer_amg_for_3D_poisson_problem(), and oomph::Biharmonic_schur_complement_Hypre_defaults::set_defaults().
| 
 | inline | 
Access function to AMG_strength.
Definition at line 1038 of file hypre_solver.h.
References oomph::HypreInterface::AMG_strength.
Referenced by oomph::Lagrange_Enforced_Flow_Preconditioner_Subsidiary_Operator_Helper::boomer_amg2v22_for_3D_momentum(), oomph::Lagrange_Enforced_Flow_Preconditioner_Subsidiary_Operator_Helper::boomer_amg_for_2D_momentum_simple_visc(), oomph::Lagrange_Enforced_Flow_Preconditioner_Subsidiary_Operator_Helper::boomer_amg_for_2D_momentum_stressdiv_visc(), oomph::Lagrange_Enforced_Flow_Preconditioner_Subsidiary_Operator_Helper::boomer_amg_for_2D_poisson_problem(), oomph::Lagrange_Enforced_Flow_Preconditioner_Subsidiary_Operator_Helper::boomer_amg_for_3D_momentum(), oomph::Lagrange_Enforced_Flow_Preconditioner_Subsidiary_Operator_Helper::boomer_amg_for_3D_poisson_problem(), oomph::Biharmonic_schur_complement_Hypre_defaults::set_defaults(), oomph::Hypre_default_settings::set_defaults_for_2D_poisson_problem(), oomph::Hypre_default_settings::set_defaults_for_3D_poisson_problem(), and oomph::Hypre_default_settings::set_defaults_for_navier_stokes_momentum_block().
| 
 | inline | 
Access function to AMG_truncation.
Definition at line 1050 of file hypre_solver.h.
References oomph::HypreInterface::AMG_truncation.
| 
 | inline | 
Function to select use of 'complex' AMG smoothers as controlled by the flag AMG_complex_smoother.
Definition at line 990 of file hypre_solver.h.
References oomph::HypreInterface::AMG_using_simple_smoothing.
| 
 | inline | 
Function to select use of 'simple' AMG smoothers as controlled by the flag AMG_simple_smoother.
Definition at line 977 of file hypre_solver.h.
References oomph::HypreInterface::AMG_using_simple_smoothing.
Referenced by oomph::Pseudo_Elastic_Preconditioner_Subsidiary_Operator_Helper::get_elastic_preconditioner_hypre(), and oomph::Hypre_default_settings::set_defaults_for_2D_poisson_problem().
| 
 | inline | 
Return function for the AMG_using_simple_smoothing_flag.
Definition at line 1002 of file hypre_solver.h.
References oomph::HypreInterface::AMG_using_simple_smoothing.
| 
 | virtual | 
Function deletes all solver data.
clean_up_memory() deletes any existing Hypre solver and Hypre matrix
Reimplemented from oomph::Preconditioner.
Definition at line 1842 of file hypre_solver.cc.
References oomph::HypreInterface::hypre_clean_up_memory().
Referenced by setup().
| 
 | inline | 
Hypre copies matrix data from oomph-lib's CRDoubleMatrix or DistributedCRDoubleMatrix into its own data structures, doubling the memory requirements for the matrix. Calling this function ensures that the matrix is not deleted.
Definition at line 1171 of file hypre_solver.h.
References Delete_matrix.
| 
 | inline | 
Disable documentation of preconditioner timings.
Definition at line 939 of file hypre_solver.h.
References Doc_time.
| 
 | inline | 
Disable euclid row scaling.
Definition at line 1110 of file hypre_solver.h.
References oomph::HypreInterface::Euclid_rowScale.
| 
 | inline | 
Disable use of Block Jacobi, so PILU will be used.
Definition at line 1124 of file hypre_solver.h.
References oomph::HypreInterface::Euclid_using_BJ.
| 
 | inline | 
Disable reporting of cumulative solve time in destructor.
Definition at line 927 of file hypre_solver.h.
References Report_my_cumulative_preconditioner_solve_time.
| 
 | inline | 
Hypre copies matrix data from oomph-lib's CRDoubleMatrix or DistributedCRDoubleMatrix into its own data structures, doubling the memory requirements for the matrix. 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 calling this function which changes the flag from false (its default) to true.
Definition at line 1162 of file hypre_solver.h.
References Delete_matrix.
| 
 | inline | 
Enable documentation of preconditioner timings.
Definition at line 933 of file hypre_solver.h.
References Doc_time.
| 
 | inline | 
Enable euclid rowScaling.
Definition at line 1104 of file hypre_solver.h.
References oomph::HypreInterface::Euclid_rowScale.
| 
 | inline | 
Enable use of Block Jacobi as opposed to PILU.
Definition at line 1117 of file hypre_solver.h.
References oomph::HypreInterface::Euclid_using_BJ.
| 
 | inline | 
Enable reporting of cumulative solve time in destructor.
Definition at line 921 of file hypre_solver.h.
References Report_my_cumulative_preconditioner_solve_time.
| 
 | inline | 
Access function to Euclid drop tolerance parameter.
Definition at line 1092 of file hypre_solver.h.
References oomph::HypreInterface::Euclid_droptol.
| 
 | inline | 
Access function to Euclid level parameter for ILU(k) factorization.
Definition at line 1098 of file hypre_solver.h.
References oomph::HypreInterface::Euclid_level.
| 
 | inline | 
Function 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 1148 of file hypre_solver.h.
References oomph::HypreInterface::Euclid_print_level.
| 
 | inline | 
Function to switch on ILU(k) factorization for Euclid (default is ILU(k) factorization)
Definition at line 1131 of file hypre_solver.h.
References oomph::HypreInterface::Euclid_using_ILUT.
| 
 | inline | 
Function to switch on ILUT factorization for Euclid (default is ILU(k) factorization)
Definition at line 1138 of file hypre_solver.h.
References oomph::HypreInterface::Euclid_using_ILUT.
| 
 | inline | 
Access function to Hypre_method flag – specified via enumeration.
Definition at line 945 of file hypre_solver.h.
References oomph::HypreInterface::Hypre_method.
Referenced by oomph::Pseudo_Elastic_Preconditioner_Subsidiary_Operator_Helper::get_elastic_preconditioner_hypre(), oomph::Biharmonic_schur_complement_Hypre_defaults::set_defaults(), and oomph::Hypre_default_settings::set_defaults_for_2D_poisson_problem().
| 
 | inline | 
Access function to Internal_preconditioner flag – specified via enumeration.
Definition at line 952 of file hypre_solver.h.
References oomph::HypreInterface::Internal_preconditioner.
| 
 | delete | 
Broken assignment operator.
| 
 | inline | 
Access function to ParaSails filter parameter.
Definition at line 1080 of file hypre_solver.h.
References oomph::HypreInterface::ParaSails_filter.
| 
 | inline | 
Access function to ParaSails nlevel parameter.
Definition at line 1068 of file hypre_solver.h.
References oomph::HypreInterface::ParaSails_nlevel.
| 
 | inline | 
Access function to ParaSails symmetry flag.
Definition at line 1062 of file hypre_solver.h.
References oomph::HypreInterface::ParaSails_symmetry.
| 
 | inline | 
Access function to ParaSails thresh parameter.
Definition at line 1074 of file hypre_solver.h.
References oomph::HypreInterface::ParaSails_thresh.
| 
 | virtual | 
Function applies solver to vector r for preconditioning. This requires a call to setup(...) first. Note: Hypre copies matrix data from oomph-lib's CRDoubleMatrix or DistributedCRDoubleMatrix into its own data structures, doubling the memory requirements for the matrix. 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 with the delete_matrix() function.
Preconditioner_solve uses a hypre solver to precondition vector r.
Implements oomph::Preconditioner.
Definition at line 1775 of file hypre_solver.cc.
References oomph::DoubleVector::built(), Context_based_cumulative_npreconditioner_solve, Context_based_cumulative_solve_time, Context_string, Cumulative_npreconditioner_solve, Cumulative_preconditioner_solve_time, oomph::DistributableLinearAlgebraObject::distribution_pt(), oomph::HypreInterface::existing_solver(), oomph::HypreInterface::hypre_solve(), My_cumulative_preconditioner_solve_time, oomph::HypreInterface::None, oomph::HypreInterface::Output_info, and oomph::TimingHelpers::timer().
| 
 | static | 
Report cumulative solve times of all instantiations of this class.
Definition at line 1648 of file hypre_solver.cc.
References Context_based_cumulative_npreconditioner_solve, Context_based_cumulative_solve_time, Context_based_nrow, Cumulative_npreconditioner_solve, Cumulative_preconditioner_solve_time, and oomph::oomph_info.
| 
 | static | 
Reset cumulative solve times.
Definition at line 1699 of file hypre_solver.cc.
References Context_based_cumulative_npreconditioner_solve, Context_based_cumulative_solve_time, Context_based_nrow, Cumulative_npreconditioner_solve, and Cumulative_preconditioner_solve_time.
| 
 | inline | 
Function to set the number of times to apply BoomerAMG.
Definition at line 964 of file hypre_solver.h.
References amg_iterations(), and oomph::HypreInterface::Max_iter.
Referenced by oomph::Lagrange_Enforced_Flow_Preconditioner_Subsidiary_Operator_Helper::boomer_amg2v22_for_3D_momentum(), oomph::Lagrange_Enforced_Flow_Preconditioner_Subsidiary_Operator_Helper::boomer_amg_for_2D_momentum_simple_visc(), oomph::Lagrange_Enforced_Flow_Preconditioner_Subsidiary_Operator_Helper::boomer_amg_for_2D_momentum_stressdiv_visc(), oomph::Lagrange_Enforced_Flow_Preconditioner_Subsidiary_Operator_Helper::boomer_amg_for_2D_poisson_problem(), oomph::Lagrange_Enforced_Flow_Preconditioner_Subsidiary_Operator_Helper::boomer_amg_for_3D_momentum(), oomph::Lagrange_Enforced_Flow_Preconditioner_Subsidiary_Operator_Helper::boomer_amg_for_3D_poisson_problem(), oomph::Pseudo_Elastic_Preconditioner_Subsidiary_Operator_Helper::get_elastic_preconditioner_hypre(), and oomph::Hypre_default_settings::set_defaults_for_2D_poisson_problem().
| 
 | virtual | 
Function to set up a preconditioner for the linear system defined by matrix_pt. This function is required when preconditioning and must be called before using the preconditioner_solve(...) function. This interface allows HyprePreconditioner to be used as a Preconditioner object for oomph-lib's own IterativeLinearSolver class. Note: matrix_pt must point to an object of type CRDoubleMatrix or DistributedCRDoubleMatrix. Note: Hypre copies matrix data from oomph-lib's CRDoubleMatrix and DistributedCRDoubleMatrix into its own data structures, doubling the memory requirements for the matrix. 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 calling the enable_delete_matrix() function.
An interface to allow HypreSolver to be used as a Preconditioner for the oomph-lib IterativeLinearSolver class. Matrix has to be of type CRDoubleMatrix or DistributedCRDoubleMatrix.
Implements oomph::Preconditioner.
Definition at line 1714 of file hypre_solver.cc.
References oomph::DistributableLinearAlgebraObject::build_distribution(), clean_up_memory(), Context_based_nrow, Context_string, oomph::HypreInterface::Delete_input_data, Delete_matrix, oomph::DistributableLinearAlgebraObject::distribution_pt(), Doc_time, oomph::HypreInterface::hypre_matrix_setup(), oomph::HypreInterface::hypre_solver_setup(), oomph::Preconditioner::matrix_pt(), oomph::DoubleMatrixBase::ncol(), oomph::CRDoubleMatrix::nnz(), oomph::DistributableLinearAlgebraObject::nrow(), oomph::CRDoubleMatrix::nrow(), oomph::DoubleMatrixBase::nrow(), oomph::DistributableLinearAlgebraObject::nrow_local(), oomph::oomph_info, and oomph::HypreInterface::Output_info.
| 
 | inline | 
Function to select BoomerAMG as the preconditioner.
Definition at line 958 of file hypre_solver.h.
References oomph::HypreInterface::BoomerAMG, and oomph::HypreInterface::Hypre_method.
| 
 | inline | 
Function to select use Euclid as the preconditioner.
Definition at line 1086 of file hypre_solver.h.
References oomph::HypreInterface::Euclid, and oomph::HypreInterface::Hypre_method.
| 
 | inline | 
Function to select ParaSails as the preconditioner.
Definition at line 1056 of file hypre_solver.h.
References oomph::HypreInterface::Hypre_method, and oomph::HypreInterface::ParaSails.
| 
 | static | 
Static unsigned that accumulates the number of preconditioner solves of all instantiations of this class, labeled by context string. Reset it manually, e.g. after every Newton solve, using reset_cumulative_solve_times().
Definition at line 904 of file hypre_solver.h.
Referenced by preconditioner_solve(), report_cumulative_solve_times(), and reset_cumulative_solve_times().
| 
 | static | 
Static double that accumulates the preconditioner solve time of all instantiations of this class, labeled by context string. Reset it manually, e.g. after every Newton solve, using reset_cumulative_solve_times().
map of static doubles that accumulates the preconditioner solve time of all instantiations of this class, labeled by context string. Reset it manually, e.g. after every Newton solve, using reset_cumulative_solve_times().
Definition at line 889 of file hypre_solver.h.
Referenced by preconditioner_solve(), report_cumulative_solve_times(), and reset_cumulative_solve_times().
| 
 | static | 
Static unsigned that stores nrow for the most recent instantiations of this class, labeled by context string. Reset it manually, e.g. after every Newton solve, using reset_cumulative_solve_times().
Definition at line 911 of file hypre_solver.h.
Referenced by report_cumulative_solve_times(), reset_cumulative_solve_times(), and setup().
| 
 | private | 
String can be used to provide context for annotation.
Definition at line 1234 of file hypre_solver.h.
Referenced by HyprePreconditioner(), preconditioner_solve(), setup(), and ~HyprePreconditioner().
| 
 | static | 
Static unsigned that accumulates the number of preconditioner solves of all instantiations of this class. Reset it manually, e.g. after every Newton solve, using reset_cumulative_solve_times().
Definition at line 896 of file hypre_solver.h.
Referenced by preconditioner_solve(), report_cumulative_solve_times(), and reset_cumulative_solve_times().
| 
 | static | 
Static double that accumulates the preconditioner solve time of all instantiations of this class. Reset it manually, e.g. after every Newton solve, using reset_cumulative_solve_times().
////////////////////////////////////////////////////////////////////
//////////////////////////////////////////////////////////////////// Static double that accumulates the preconditioner solve time of all instantiations of this class. Reset it manually, e.g. after every Newton solve, using reset_cumulative_solve_times().
Definition at line 882 of file hypre_solver.h.
Referenced by preconditioner_solve(), report_cumulative_solve_times(), and reset_cumulative_solve_times().
| 
 | private | 
Hypre copies matrix data from oomph-lib's CRDoubleMatrix or DistributedCRDoubleMatrix into its own data structures, doubling the memory requirements for the matrix. 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 1218 of file hypre_solver.h.
Referenced by disable_delete_matrix(), enable_delete_matrix(), HyprePreconditioner(), and setup().
| 
 | private | 
Definition at line 1221 of file hypre_solver.h.
Referenced by disable_doc_time(), enable_doc_time(), HyprePreconditioner(), and setup().
| 
 | private | 
Private double that accumulates the preconditioner solve time of thi instantiation of this class. Is reported in destructor if Report_my_cumulative_preconditioner_solve_time is set to true.
Definition at line 1227 of file hypre_solver.h.
Referenced by HyprePreconditioner(), preconditioner_solve(), and ~HyprePreconditioner().
| 
 | private | 
Bool to request report of My_cumulative_preconditioner_solve_time in destructor.
Definition at line 1231 of file hypre_solver.h.
Referenced by disable_report_my_cumulative_preconditioner_solve_time(), enable_report_my_cumulative_preconditioner_solve_time(), HyprePreconditioner(), and ~HyprePreconditioner().