Public Member Functions | Private Attributes | List of all members
oomph::HypreSolver Class Reference

An LinearSolver class using the suite of Hypre solvers to allow. More...

#include <hypre_solver.h>

+ Inheritance diagram for oomph::HypreSolver:

Public Member Functions

 HypreSolver ()
 Constructor.
 
 ~HypreSolver ()
 Empty destructor.
 
 HypreSolver (const HypreSolver &)=delete
 Broken copy constructor.
 
void operator= (const HypreSolver &)=delete
 Broken assignment operator.
 
void disable_resolve ()
 Disable resolve function (overloads the LinearSolver disable_resolve function).
 
unsignedmax_iter ()
 Access function to Max_iter.
 
doubletolerance ()
 Access function to Tolerance.
 
unsignedhypre_method ()
 Access function to Hypre_method flag – specified via enumeration.
 
unsignedinternal_preconditioner ()
 Access function to Internal_preconditioner flag – specified via enumeration.
 
void amg_using_simple_smoothing ()
 Function to select use of 'simple' AMG smoothers as controlled by AMG_simple_smoother flag.
 
unsignedamg_simple_smoother ()
 Access function to AMG_simple_smoother flag.
 
void amg_using_complex_smoothing ()
 Function to select use of 'complex' AMG smoothers as controlled by AMG_complex_smoother flag.
 
unsignedamg_complex_smoother ()
 Access function to AMG_complex_smoother flag.
 
unsignedamg_print_level ()
 Access function to AMG_print_level.
 
unsignedamg_smoother_iterations ()
 Access function to AMG_smoother_iterations.
 
unsignedamg_coarsening ()
 Access function to AMG_coarsening flag.
 
unsignedamg_max_levels ()
 Access function to AMG_max_levels.
 
doubleamg_damping ()
 Access function to AMG_damping parameter.
 
doubleamg_strength ()
 Access function to AMG_strength.
 
doubleamg_max_row_sum ()
 Access function to AMG_max_row_sum.
 
doubleamg_truncation ()
 Access function to AMG_truncation.
 
intparasails_symmetry ()
 Access function to ParaSails symmetry flag.
 
intparasails_nlevel ()
 Access function to ParaSails nlevel parameter.
 
doubleparasails_thresh ()
 Access function to ParaSails thresh parameter.
 
doubleparasails_filter ()
 Access function to ParaSails filter parameter.
 
doubleeuclid_droptol ()
 Access function to Euclid drop tolerance parameter.
 
inteuclid_level ()
 Access function to Euclid level parameter for ILU(k) factorization.
 
void enable_euclid_rowScale ()
 Enable euclid rowScaling.
 
void disable_euclid_rowScale ()
 Disable euclid row scaling.
 
void enable_euclid_using_BJ ()
 Enable use of Block Jacobi as opposed to PILU.
 
void disable_euclid_using_BJ ()
 Disable use of Block Jacobi, so PILU will be used.
 
void euclid_using_ILUK ()
 Function to switch on ILU(k) factorization for Euclid (default is ILU(k) factorization)
 
void euclid_using_ILUT ()
 Function to switch on ILUT factorization for Euclid (default is ILU(k) factorization)
 
unsignedeuclid_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.
 
unsignedkrylov_print_level ()
 Access function to Krylov_print_level.
 
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.
 
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.
 
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 solved. This function deletes any existing internal data and then generates a new Hypre solver.
 
void solve (DoubleMatrixBase *const &matrix_pt, const DoubleVector &rhs, DoubleVector &solution)
 Function to solve the linear system defined by matrix_pt and rhs. This function will delete any existing internal data and generate a new Hypre solver. Note: The matrix has to be of type CRDoubleMatrix or Distributed CRDoubleMatrix. 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 by calling the enable_delete_matrix() function.
 
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 hand side vector. This function must be used after a call to solve(...) with enable_resolve set to true.
 
void clean_up_memory ()
 Function deletes all solver data.
 
- Public Member Functions inherited from oomph::LinearSolver
 LinearSolver ()
 Empty constructor, initialise the member data.
 
 LinearSolver (const LinearSolver &dummy)=delete
 Broken copy constructor.
 
void operator= (const LinearSolver &)=delete
 Broken assignment operator.
 
virtual ~LinearSolver ()
 Empty virtual destructor.
 
void enable_doc_time ()
 Enable documentation of solve times.
 
void disable_doc_time ()
 Disable documentation of solve times.
 
bool is_doc_time_enabled () const
 Is documentation of solve times enabled?
 
bool is_resolve_enabled () const
 Boolean flag indicating if resolves are enabled.
 
virtual void enable_resolve ()
 Enable resolve (i.e. store matrix and/or LU decomposition, say) Virtual so it can be overloaded to perform additional tasks.
 
virtual void solve (DoubleMatrixBase *const &matrix_pt, const Vector< double > &rhs, Vector< double > &result)
 Linear-algebra-type solver: Takes pointer to a matrix and rhs vector and returns the solution of the linear system.
 
virtual void solve_transpose (Problem *const &problem_pt, DoubleVector &result)
 Solver: Takes pointer to problem and returns the results vector which contains the solution of the linear system defined by the problem's fully assembled Jacobian and residual vector (broken virtual).
 
virtual void solve_transpose (DoubleMatrixBase *const &matrix_pt, const DoubleVector &rhs, DoubleVector &result)
 Linear-algebra-type solver: Takes pointer to a matrix and rhs vector and returns the solution of the linear system.
 
virtual void solve_transpose (DoubleMatrixBase *const &matrix_pt, const Vector< double > &rhs, Vector< double > &result)
 Linear-algebra-type solver: Takes pointer to a matrix and rhs vector and returns the solution of the linear system.
 
virtual void resolve_transpose (const DoubleVector &rhs, DoubleVector &result)
 Solver: Resolve the system defined by the last assembled jacobian and the rhs vector. Solution is returned in the vector result. (broken virtual)
 
virtual double jacobian_setup_time () const
 returns the time taken to assemble the Jacobian matrix and residual vector (needs to be overloaded for each solver)
 
virtual double linear_solver_solution_time () const
 return the time taken to solve the linear system (needs to be overloaded for each linear solver)
 
virtual void enable_computation_of_gradient ()
 function to enable the computation of the gradient required for the globally convergent Newton method
 
void disable_computation_of_gradient ()
 function to disable the computation of the gradient required for the globally convergent Newton method
 
void reset_gradient ()
 function to reset the size of the gradient before each Newton solve
 
void get_gradient (DoubleVector &gradient)
 function to access the gradient, provided it has been computed
 
- Public Member Functions inherited from oomph::DistributableLinearAlgebraObject
 DistributableLinearAlgebraObject ()
 Default constructor - create a distribution.
 
 DistributableLinearAlgebraObject (const DistributableLinearAlgebraObject &matrix)=delete
 Broken copy constructor.
 
void operator= (const DistributableLinearAlgebraObject &)=delete
 Broken assignment operator.
 
virtual ~DistributableLinearAlgebraObject ()
 Destructor.
 
LinearAlgebraDistributiondistribution_pt () const
 access to the LinearAlgebraDistribution
 
unsigned nrow () const
 access function to the number of global rows.
 
unsigned nrow_local () const
 access function for the num of local rows on this processor.
 
unsigned nrow_local (const unsigned &p) const
 access function for the num of local rows on this processor.
 
unsigned first_row () const
 access function for the first row on this processor
 
unsigned first_row (const unsigned &p) const
 access function for the first row on this processor
 
bool distributed () const
 distribution is serial or distributed
 
bool distribution_built () const
 if the communicator_pt is null then the distribution is not setup then false is returned, otherwise return true
 
void build_distribution (const LinearAlgebraDistribution *const dist_pt)
 setup the distribution of this distributable linear algebra object
 
void build_distribution (const LinearAlgebraDistribution &dist)
 setup the distribution of this distributable linear algebra object
 
- Public Member Functions inherited from oomph::HypreInterface
 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.
 

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.
 

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...
 
- Protected Member Functions inherited from oomph::DistributableLinearAlgebraObject
void clear_distribution ()
 clear the distribution of this distributable linear algebra object
 
- Protected Member Functions inherited from oomph::HypreInterface
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 inherited from oomph::LinearSolver
bool Enable_resolve
 Boolean that indicates whether the matrix (or its factors, in the case of direct solver) should be stored so that the resolve function can be used.
 
bool Doc_time
 Boolean flag that indicates whether the time taken.
 
bool Compute_gradient
 flag that indicates whether the gradient required for the globally convergent Newton method should be computed or not
 
bool Gradient_has_been_computed
 flag that indicates whether the gradient was computed or not
 
DoubleVector Gradient_for_glob_conv_newton_solve
 DoubleVector storing the gradient for the globally convergent Newton method.
 
- Protected Attributes inherited from oomph::HypreInterface
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.
 

Detailed Description

An LinearSolver class using the suite of Hypre solvers to allow.

BoomerAMG (AMG), CG, GMRES or BiCGStab

to be used for LinearSolver::solve(...) or LinearSolver::resolve(...).

The Krylov subspace solvers (CG, GMRES and BiCGStab) may be preconditioned using:

    BoomerAMG, Euclid or ParaSails

By default GMRES without preconditioning is used.

Definition at line 702 of file hypre_solver.h.

Constructor & Destructor Documentation

◆ HypreSolver() [1/2]

oomph::HypreSolver::HypreSolver ( )
inline

Constructor.

Definition at line 706 of file hypre_solver.h.

References Delete_matrix, and oomph::LinearSolver::Doc_time.

◆ ~HypreSolver()

oomph::HypreSolver::~HypreSolver ( )
inline

Empty destructor.

Definition at line 723 of file hypre_solver.h.

◆ HypreSolver() [2/2]

oomph::HypreSolver::HypreSolver ( const HypreSolver )
delete

Broken copy constructor.

Member Function Documentation

◆ amg_coarsening()

unsigned & oomph::HypreSolver::amg_coarsening ( )
inline

Access function to AMG_coarsening flag.

Definition at line 803 of file hypre_solver.h.

References oomph::HypreInterface::AMG_coarsening.

◆ amg_complex_smoother()

unsigned & oomph::HypreSolver::amg_complex_smoother ( )
inline

Access function to AMG_complex_smoother flag.

Definition at line 785 of file hypre_solver.h.

References oomph::HypreInterface::AMG_complex_smoother.

◆ amg_damping()

double & oomph::HypreSolver::amg_damping ( )
inline

Access function to AMG_damping parameter.

Definition at line 815 of file hypre_solver.h.

References oomph::HypreInterface::AMG_damping.

◆ amg_max_levels()

unsigned & oomph::HypreSolver::amg_max_levels ( )
inline

Access function to AMG_max_levels.

Definition at line 809 of file hypre_solver.h.

References oomph::HypreInterface::AMG_max_levels.

◆ amg_max_row_sum()

double & oomph::HypreSolver::amg_max_row_sum ( )
inline

Access function to AMG_max_row_sum.

Definition at line 827 of file hypre_solver.h.

References oomph::HypreInterface::AMG_max_row_sum.

◆ amg_print_level()

unsigned & oomph::HypreSolver::amg_print_level ( )
inline

Access function to AMG_print_level.

Definition at line 791 of file hypre_solver.h.

References oomph::HypreInterface::AMG_print_level.

◆ amg_simple_smoother()

unsigned & oomph::HypreSolver::amg_simple_smoother ( )
inline

Access function to AMG_simple_smoother flag.

Definition at line 772 of file hypre_solver.h.

References oomph::HypreInterface::AMG_simple_smoother.

◆ amg_smoother_iterations()

unsigned & oomph::HypreSolver::amg_smoother_iterations ( )
inline

Access function to AMG_smoother_iterations.

Definition at line 797 of file hypre_solver.h.

References oomph::HypreInterface::AMG_smoother_iterations.

◆ amg_strength()

double & oomph::HypreSolver::amg_strength ( )
inline

Access function to AMG_strength.

Definition at line 821 of file hypre_solver.h.

References oomph::HypreInterface::AMG_strength.

◆ amg_truncation()

double & oomph::HypreSolver::amg_truncation ( )
inline

Access function to AMG_truncation.

Definition at line 833 of file hypre_solver.h.

References oomph::HypreInterface::AMG_truncation.

◆ amg_using_complex_smoothing()

void oomph::HypreSolver::amg_using_complex_smoothing ( )
inline

Function to select use of 'complex' AMG smoothers as controlled by AMG_complex_smoother flag.

Definition at line 779 of file hypre_solver.h.

References oomph::HypreInterface::AMG_using_simple_smoothing.

◆ amg_using_simple_smoothing()

void oomph::HypreSolver::amg_using_simple_smoothing ( )
inline

Function to select use of 'simple' AMG smoothers as controlled by AMG_simple_smoother flag.

Definition at line 766 of file hypre_solver.h.

References oomph::HypreInterface::AMG_using_simple_smoothing.

◆ clean_up_memory()

void oomph::HypreSolver::clean_up_memory ( )
virtual

Function deletes all solver data.

clean_up_memory() deletes any existing solver data

Reimplemented from oomph::LinearSolver.

Definition at line 1649 of file hypre_solver.cc.

References oomph::HypreInterface::hypre_clean_up_memory().

Referenced by disable_resolve(), solve(), and solve().

◆ disable_delete_matrix()

void oomph::HypreSolver::disable_delete_matrix ( )
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 949 of file hypre_solver.h.

References Delete_matrix.

◆ disable_euclid_rowScale()

void oomph::HypreSolver::disable_euclid_rowScale ( )
inline

Disable euclid row scaling.

Definition at line 882 of file hypre_solver.h.

References oomph::HypreInterface::Euclid_rowScale.

◆ disable_euclid_using_BJ()

void oomph::HypreSolver::disable_euclid_using_BJ ( )
inline

Disable use of Block Jacobi, so PILU will be used.

Definition at line 896 of file hypre_solver.h.

References oomph::HypreInterface::Euclid_using_BJ.

◆ disable_resolve()

void oomph::HypreSolver::disable_resolve ( )
inlinevirtual

Disable resolve function (overloads the LinearSolver disable_resolve function).

Reimplemented from oomph::LinearSolver.

Definition at line 733 of file hypre_solver.h.

References clean_up_memory(), and oomph::LinearSolver::Enable_resolve.

◆ enable_delete_matrix()

void oomph::HypreSolver::enable_delete_matrix ( )
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 940 of file hypre_solver.h.

References Delete_matrix.

◆ enable_euclid_rowScale()

void oomph::HypreSolver::enable_euclid_rowScale ( )
inline

Enable euclid rowScaling.

Definition at line 876 of file hypre_solver.h.

References oomph::HypreInterface::Euclid_rowScale.

◆ enable_euclid_using_BJ()

void oomph::HypreSolver::enable_euclid_using_BJ ( )
inline

Enable use of Block Jacobi as opposed to PILU.

Definition at line 889 of file hypre_solver.h.

References oomph::HypreInterface::Euclid_using_BJ.

◆ euclid_droptol()

double & oomph::HypreSolver::euclid_droptol ( )
inline

Access function to Euclid drop tolerance parameter.

Definition at line 863 of file hypre_solver.h.

References oomph::HypreInterface::Euclid_droptol.

◆ euclid_level()

int & oomph::HypreSolver::euclid_level ( )
inline

Access function to Euclid level parameter for ILU(k) factorization.

Definition at line 870 of file hypre_solver.h.

References oomph::HypreInterface::Euclid_level.

◆ euclid_print_level()

unsigned & oomph::HypreSolver::euclid_print_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 920 of file hypre_solver.h.

References oomph::HypreInterface::Euclid_print_level.

◆ euclid_using_ILUK()

void oomph::HypreSolver::euclid_using_ILUK ( )
inline

Function to switch on ILU(k) factorization for Euclid (default is ILU(k) factorization)

Definition at line 903 of file hypre_solver.h.

References oomph::HypreInterface::Euclid_using_ILUT.

◆ euclid_using_ILUT()

void oomph::HypreSolver::euclid_using_ILUT ( )
inline

Function to switch on ILUT factorization for Euclid (default is ILU(k) factorization)

Definition at line 910 of file hypre_solver.h.

References oomph::HypreInterface::Euclid_using_ILUT.

◆ hypre_method()

unsigned & oomph::HypreSolver::hypre_method ( )
inline

Access function to Hypre_method flag – specified via enumeration.

Definition at line 752 of file hypre_solver.h.

References oomph::HypreInterface::Hypre_method.

◆ internal_preconditioner()

unsigned & oomph::HypreSolver::internal_preconditioner ( )
inline

Access function to Internal_preconditioner flag – specified via enumeration.

Definition at line 759 of file hypre_solver.h.

References oomph::HypreInterface::Internal_preconditioner.

◆ krylov_print_level()

unsigned & oomph::HypreSolver::krylov_print_level ( )
inline

Access function to Krylov_print_level.

Definition at line 926 of file hypre_solver.h.

References oomph::HypreInterface::Krylov_print_level.

◆ max_iter()

unsigned & oomph::HypreSolver::max_iter ( )
inline

Access function to Max_iter.

Definition at line 740 of file hypre_solver.h.

References oomph::HypreInterface::Max_iter.

◆ operator=()

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

Broken assignment operator.

◆ parasails_filter()

double & oomph::HypreSolver::parasails_filter ( )
inline

Access function to ParaSails filter parameter.

Definition at line 857 of file hypre_solver.h.

References oomph::HypreInterface::ParaSails_filter.

◆ parasails_nlevel()

int & oomph::HypreSolver::parasails_nlevel ( )
inline

Access function to ParaSails nlevel parameter.

Definition at line 845 of file hypre_solver.h.

References oomph::HypreInterface::ParaSails_nlevel.

◆ parasails_symmetry()

int & oomph::HypreSolver::parasails_symmetry ( )
inline

Access function to ParaSails symmetry flag.

Definition at line 839 of file hypre_solver.h.

References oomph::HypreInterface::ParaSails_symmetry.

◆ parasails_thresh()

double & oomph::HypreSolver::parasails_thresh ( )
inline

Access function to ParaSails thresh parameter.

Definition at line 851 of file hypre_solver.h.

References oomph::HypreInterface::ParaSails_thresh.

◆ resolve()

void oomph::HypreSolver::resolve ( const DoubleVector rhs,
DoubleVector solution 
)
virtual

Function to resolve a linear system using the existing solver data, allowing a solve with a new right hand side vector. This function must be used after a call to solve(...) with enable_resolve set to true.

Resolve performs a linear solve using current solver data (if such data exists).

Reimplemented from oomph::LinearSolver.

Definition at line 1595 of file hypre_solver.cc.

References oomph::DistributableLinearAlgebraObject::distribution_pt(), oomph::LinearSolver::Doc_time, oomph::HypreInterface::existing_solver(), oomph::HypreInterface::hypre_solve(), oomph::HypreInterface::None, and oomph::HypreInterface::Output_info.

◆ solve() [1/2]

void oomph::HypreSolver::solve ( DoubleMatrixBase *const matrix_pt,
const DoubleVector rhs,
DoubleVector solution 
)
virtual

Function to solve the linear system defined by matrix_pt and rhs. This function will delete any existing internal data and generate a new Hypre solver. Note: The matrix has to be of type CRDoubleMatrix or Distributed CRDoubleMatrix. 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 by calling the enable_delete_matrix() function.

Uses HypreInterface::hypre_solve(...) to solve the linear system for a CRDoubleMatrix or a DistributedCRDoubleMatrix. In the latter case, the rhs needs to be distributed too, i.e. the length of the rhs vector must be equal to the number of rows stored locally. Note: the returned solution vector is never distributed, i.e. all processors hold all values.

Reimplemented from oomph::LinearSolver.

Definition at line 1499 of file hypre_solver.cc.

References oomph::DistributableLinearAlgebraObject::build_distribution(), clean_up_memory(), oomph::HypreInterface::Delete_input_data, Delete_matrix, oomph::DistributableLinearAlgebraObject::distribution_pt(), oomph::LinearSolver::Doc_time, oomph::LinearSolver::Enable_resolve, oomph::HypreInterface::hypre_matrix_setup(), oomph::HypreInterface::hypre_solve(), oomph::HypreInterface::hypre_solver_setup(), oomph::DoubleMatrixBase::ncol(), oomph::DoubleMatrixBase::nrow(), and oomph::HypreInterface::Output_info.

◆ solve() [2/2]

void oomph::HypreSolver::solve ( Problem *const problem_pt,
DoubleVector solution 
)
virtual

Function which uses problem_pt's get_jacobian(...) function to generate a linear system which is then solved. This function deletes any existing internal data and then generates a new Hypre solver.

Problem-based solve function to generate the Jacobian matrix and residual vector and use HypreInterface::hypre_solver_setup(...) and HypreInterface::hypre_solve(...) to solve the linear system. This function will delete any existing data. Note: the returned solution vector is NOT distributed, i.e. all processors hold all values because this is what the Newton solver requires.

Implements oomph::LinearSolver.

Definition at line 1446 of file hypre_solver.cc.

References clean_up_memory(), oomph::HypreInterface::Delete_input_data, oomph::LinearSolver::Doc_time, oomph::LinearSolver::Enable_resolve, oomph::Problem::get_jacobian(), oomph::HypreInterface::hypre_matrix_setup(), oomph::HypreInterface::hypre_solve(), oomph::HypreInterface::hypre_solver_setup(), oomph::oomph_info, oomph::HypreInterface::Output_info, and oomph::TimingHelpers::timer().

◆ tolerance()

double & oomph::HypreSolver::tolerance ( )
inline

Access function to Tolerance.

Definition at line 746 of file hypre_solver.h.

References oomph::HypreInterface::Tolerance.

Member Data Documentation

◆ Delete_matrix

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

Referenced by disable_delete_matrix(), enable_delete_matrix(), HypreSolver(), and solve().


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