oomph::LinearSolver Class Referenceabstract

#include <linear_solver.h>

+ Inheritance diagram for oomph::LinearSolver:

Public Member Functions

 LinearSolver ()
 Empty constructor, initialise the member data. More...
 
 LinearSolver (const LinearSolver &dummy)=delete
 Broken copy constructor. More...
 
void operator= (const LinearSolver &)=delete
 Broken assignment operator. More...
 
virtual ~LinearSolver ()
 Empty virtual destructor. More...
 
void enable_doc_time ()
 Enable documentation of solve times. More...
 
void disable_doc_time ()
 Disable documentation of solve times. More...
 
bool is_doc_time_enabled () const
 Is documentation of solve times enabled? More...
 
bool is_resolve_enabled () const
 Boolean flag indicating if resolves are enabled. More...
 
virtual void enable_resolve ()
 
virtual void disable_resolve ()
 
virtual void solve (Problem *const &problem_pt, DoubleVector &result)=0
 
virtual void solve (DoubleMatrixBase *const &matrix_pt, const DoubleVector &rhs, DoubleVector &result)
 
virtual void solve (DoubleMatrixBase *const &matrix_pt, const Vector< double > &rhs, Vector< double > &result)
 
virtual void solve_transpose (Problem *const &problem_pt, DoubleVector &result)
 
virtual void solve_transpose (DoubleMatrixBase *const &matrix_pt, const DoubleVector &rhs, DoubleVector &result)
 
virtual void solve_transpose (DoubleMatrixBase *const &matrix_pt, const Vector< double > &rhs, Vector< double > &result)
 
virtual void resolve (const DoubleVector &rhs, DoubleVector &result)
 
virtual void resolve_transpose (const DoubleVector &rhs, DoubleVector &result)
 
virtual void clean_up_memory ()
 
virtual double jacobian_setup_time () const
 
virtual double linear_solver_solution_time () const
 
virtual void enable_computation_of_gradient ()
 
void disable_computation_of_gradient ()
 
void reset_gradient ()
 
void get_gradient (DoubleVector &gradient)
 function to access the gradient, provided it has been computed 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...
 
LinearAlgebraDistributiondistribution_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
 
void build_distribution (const LinearAlgebraDistribution *const dist_pt)
 
void build_distribution (const LinearAlgebraDistribution &dist)
 

Protected Attributes

bool Enable_resolve
 
bool Doc_time
 Boolean flag that indicates whether the time taken. More...
 
bool Compute_gradient
 
bool Gradient_has_been_computed
 flag that indicates whether the gradient was computed or not More...
 
DoubleVector Gradient_for_glob_conv_newton_solve
 

Additional Inherited Members

- Protected Member Functions inherited from oomph::DistributableLinearAlgebraObject
void clear_distribution ()
 

Detailed Description

Base class for all linear solvers. This merely defines standard interfaces for linear solvers, so that different solvers can be used in a clean and transparent manner. Note that LinearSolvers are primarily used to solve the linear systems arising in oomph-lib's Newton iteration. Their primary solve function therefore takes a pointer to the associated problem, construct its Jacobian matrix and residual vector, and return the solution of the linear system formed by the Jacobian and the residual vector. We also provide broken virtual interfaces to a linear-algebra-type solve function in which the matrix and the rhs can be specified, but this are not guaranteed to implemented for all linear solvers (e.g. for frontal solvers).

Constructor & Destructor Documentation

◆ LinearSolver() [1/2]

oomph::LinearSolver::LinearSolver ( )
inline

Empty constructor, initialise the member data.

93  : Enable_resolve(false),
94  Doc_time(true),
95  Compute_gradient(false),
97  {
98  }
bool Doc_time
Boolean flag that indicates whether the time taken.
Definition: linear_solver.h:77
bool Enable_resolve
Definition: linear_solver.h:73
bool Gradient_has_been_computed
flag that indicates whether the gradient was computed or not
Definition: linear_solver.h:84
bool Compute_gradient
Definition: linear_solver.h:81

◆ LinearSolver() [2/2]

oomph::LinearSolver::LinearSolver ( const LinearSolver dummy)
delete

Broken copy constructor.

◆ ~LinearSolver()

virtual oomph::LinearSolver::~LinearSolver ( )
inlinevirtual

Empty virtual destructor.

107 {}

Member Function Documentation

◆ clean_up_memory()

◆ disable_computation_of_gradient()

void oomph::LinearSolver::disable_computation_of_gradient ( )
inline

function to disable the computation of the gradient required for the globally convergent Newton method

294  {
295  Compute_gradient = false;
296  }

References Compute_gradient.

◆ disable_doc_time()

◆ disable_resolve()

virtual void oomph::LinearSolver::disable_resolve ( )
inlinevirtual

Disable resolve (i.e. store matrix and/or LU decomposition, say) This function simply resets an internal flag. It's virtual so it can be overloaded to perform additional tasks such as cleaning up memory that is only required for the resolve.

Reimplemented in oomph::HelmholtzGMRESMG< MATRIX >, oomph::ComplexGMRES< MATRIX >, oomph::TrilinosAztecOOSolver, oomph::MumpsSolver, oomph::SuperLUSolver, oomph::AugmentedProblemGMRES, oomph::GMRES< MATRIX >, oomph::GS< CRDoubleMatrix >, oomph::GS< MATRIX >, oomph::BiCGStab< MATRIX >, oomph::CG< MATRIX >, oomph::HypreSolver, and oomph::HSL_MA42.

145  {
146  Enable_resolve = false;
147  }

References Enable_resolve.

Referenced by oomph::Problem::calculate_continuation_derivatives(), oomph::HSL_MA42::disable_resolve(), oomph::CG< MATRIX >::disable_resolve(), oomph::BiCGStab< MATRIX >::disable_resolve(), oomph::GS< MATRIX >::disable_resolve(), oomph::GS< CRDoubleMatrix >::disable_resolve(), oomph::GMRES< MATRIX >::disable_resolve(), oomph::AugmentedProblemGMRES::disable_resolve(), oomph::SuperLUSolver::disable_resolve(), oomph::MumpsSolver::disable_resolve(), oomph::ComplexGMRES< MATRIX >::disable_resolve(), oomph::HelmholtzGMRESMG< MATRIX >::disable_resolve(), oomph::FoldHandler::FoldHandler(), oomph::HopfHandler::HopfHandler(), oomph::Problem::newton_solve_continuation(), oomph::AugmentedBlockFoldLinearSolver::resolve(), oomph::BlockPitchForkLinearSolver::resolve(), oomph::AugmentedBlockPitchForkLinearSolver::resolve(), run_it(), oomph::AugmentedBlockFoldLinearSolver::solve(), oomph::BlockPitchForkLinearSolver::solve(), oomph::AugmentedBlockPitchForkLinearSolver::solve(), oomph::BlockHopfLinearSolver::solve(), oomph::ARPACK::solve_eigenproblem(), and oomph::BlockHopfLinearSolver::solve_for_two_rhs().

◆ enable_computation_of_gradient()

virtual void oomph::LinearSolver::enable_computation_of_gradient ( )
inlinevirtual

function to enable the computation of the gradient required for the globally convergent Newton method

Reimplemented in oomph::SuperLUSolver, and oomph::GMRES< MATRIX >.

283  {
284  throw OomphLibError(
285  "enable_computation_of_gradient() not implemented for "
286  "this linear solver",
289  }
#define OOMPH_EXCEPTION_LOCATION
Definition: oomph_definitions.h:61
#define OOMPH_CURRENT_FUNCTION
Definition: oomph_definitions.h:86

References OOMPH_CURRENT_FUNCTION, and OOMPH_EXCEPTION_LOCATION.

Referenced by oomph::Problem::newton_solve().

◆ enable_doc_time()

void oomph::LinearSolver::enable_doc_time ( )
inline

Enable documentation of solve times.

111  {
112  Doc_time = true;
113  }

References Doc_time.

Referenced by oomph::NewMumpsPreconditioner::enable_doc_time(), main(), and oomph::ProjectionProblem< PROJECTABLE_ELEMENT >::project().

◆ enable_resolve()

◆ get_gradient()

void oomph::LinearSolver::get_gradient ( DoubleVector gradient)
inline

function to access the gradient, provided it has been computed

307  {
308 #ifdef PARANOID
310  {
311 #endif
313 #ifdef PARANOID
314  }
315  else
316  {
317  throw OomphLibError(
318  "The gradient has not been computed for this linear solver!",
321  }
322 #endif
323  }
DoubleVector Gradient_for_glob_conv_newton_solve
Definition: linear_solver.h:88

References Gradient_for_glob_conv_newton_solve, Gradient_has_been_computed, OOMPH_CURRENT_FUNCTION, and OOMPH_EXCEPTION_LOCATION.

Referenced by oomph::Problem::newton_solve().

◆ is_doc_time_enabled()

bool oomph::LinearSolver::is_doc_time_enabled ( ) const
inline

Is documentation of solve times enabled?

123  {
124  return Doc_time;
125  }

References Doc_time.

Referenced by oomph::ProjectionProblem< PROJECTABLE_ELEMENT >::project().

◆ is_resolve_enabled()

bool oomph::LinearSolver::is_resolve_enabled ( ) const
inline

◆ jacobian_setup_time()

virtual double oomph::LinearSolver::jacobian_setup_time ( ) const
inlinevirtual

returns the time taken to assemble the Jacobian matrix and residual vector (needs to be overloaded for each solver)

Reimplemented in oomph::SuperLUSolver, oomph::DenseLU, and oomph::IterativeLinearSolver.

261  {
262  throw OomphLibError(
263  "jacobian_setup_time has not been implemented for this linear solver",
266  return 0;
267  }

References OOMPH_CURRENT_FUNCTION, and OOMPH_EXCEPTION_LOCATION.

◆ linear_solver_solution_time()

virtual double oomph::LinearSolver::linear_solver_solution_time ( ) const
inlinevirtual

return the time taken to solve the linear system (needs to be overloaded for each linear solver)

Reimplemented in oomph::SuperLUSolver, oomph::DenseLU, and oomph::IterativeLinearSolver.

272  {
273  throw OomphLibError(
274  "linear_solver_solution_time() not implemented for this linear solver",
277  return 0;
278  }

References OOMPH_CURRENT_FUNCTION, and OOMPH_EXCEPTION_LOCATION.

◆ operator=()

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

Broken assignment operator.

◆ reset_gradient()

void oomph::LinearSolver::reset_gradient ( )
inline

function to reset the size of the gradient before each Newton solve

301  {
303  }
void clear()
wipes the DoubleVector
Definition: double_vector.h:142

References oomph::DoubleVector::clear(), and Gradient_for_glob_conv_newton_solve.

Referenced by oomph::Problem::newton_solve().

◆ resolve()

virtual void oomph::LinearSolver::resolve ( const DoubleVector rhs,
DoubleVector result 
)
inlinevirtual

Resolve the system defined by the last assembled jacobian and the rhs vector. Solution is returned in the vector result. (broken virtual)

Reimplemented in oomph::TrilinosAztecOOSolver, oomph::HypreSolver, oomph::HelmholtzGMRESMG< MATRIX >, oomph::MumpsSolver, oomph::SuperLUSolver, oomph::AugmentedProblemGMRES, oomph::GMRES< MATRIX >, oomph::DampedJacobi< MATRIX >, oomph::GS< CRDoubleMatrix >, oomph::GS< MATRIX >, oomph::BiCGStab< MATRIX >, oomph::CG< MATRIX >, oomph::HSL_MA42, oomph::BlockHopfLinearSolver, oomph::AugmentedBlockPitchForkLinearSolver, oomph::BlockPitchForkLinearSolver, and oomph::AugmentedBlockFoldLinearSolver.

226  {
227  throw OomphLibError(
228  "Resolve function not implemented for this linear solver",
231  }

References OOMPH_CURRENT_FUNCTION, and OOMPH_EXCEPTION_LOCATION.

Referenced by oomph::ProblemBasedShiftInvertOperator::apply(), oomph::Problem::calculate_continuation_derivatives(), oomph::FoldHandler::FoldHandler(), oomph::Problem::get_inverse_mass_matrix_times_residuals(), oomph::HopfHandler::HopfHandler(), oomph::Problem::newton_solve(), oomph::Problem::newton_solve_continuation(), oomph::AugmentedBlockFoldLinearSolver::resolve(), oomph::BlockPitchForkLinearSolver::resolve(), oomph::AugmentedBlockPitchForkLinearSolver::resolve(), run_it(), oomph::AugmentedBlockFoldLinearSolver::solve(), oomph::BlockPitchForkLinearSolver::solve(), oomph::AugmentedBlockPitchForkLinearSolver::solve(), oomph::BlockHopfLinearSolver::solve(), oomph::ARPACK::solve_eigenproblem(), and oomph::BlockHopfLinearSolver::solve_for_two_rhs().

◆ resolve_transpose()

virtual void oomph::LinearSolver::resolve_transpose ( const DoubleVector rhs,
DoubleVector result 
)
inlinevirtual

Solver: Resolve the system defined by the last assembled jacobian and the rhs vector. Solution is returned in the vector result. (broken virtual)

Reimplemented in oomph::SuperLUSolver.

238  {
239  // Create an output stream
240  std::ostringstream error_message_stream;
241 
242  // Create the error message
243  error_message_stream
244  << "The function to resolve the transposed system has "
245  << "not yet been\nimplemented for this solver." << std::endl;
246 
247  // Throw the error message
248  throw OomphLibError(error_message_stream.str(),
251  } // End of resolve_transpose

References OOMPH_CURRENT_FUNCTION, and OOMPH_EXCEPTION_LOCATION.

◆ solve() [1/3]

virtual void oomph::LinearSolver::solve ( DoubleMatrixBase *const &  matrix_pt,
const DoubleVector rhs,
DoubleVector result 
)
inlinevirtual

◆ solve() [2/3]

virtual void oomph::LinearSolver::solve ( DoubleMatrixBase *const &  matrix_pt,
const Vector< double > &  rhs,
Vector< double > &  result 
)
inlinevirtual

Linear-algebra-type solver: Takes pointer to a matrix and rhs vector and returns the solution of the linear system.

Reimplemented in oomph::HelmholtzGMRESMG< MATRIX >, oomph::ComplexGMRES< MATRIX >, oomph::FD_LU, oomph::DenseLU, oomph::AugmentedProblemGMRES, oomph::GMRES< MATRIX >, oomph::GS< CRDoubleMatrix >, oomph::GS< MATRIX >, oomph::BiCGStab< MATRIX >, oomph::HSL_MA42, oomph::BlockHopfLinearSolver, oomph::AugmentedBlockPitchForkLinearSolver, oomph::BlockPitchForkLinearSolver, and oomph::AugmentedBlockFoldLinearSolver.

171  {
172  throw OomphLibError(
173  "Vector<double> based solve function not implemented for this solver",
176  }

References OOMPH_CURRENT_FUNCTION, and OOMPH_EXCEPTION_LOCATION.

◆ solve() [3/3]

virtual void oomph::LinearSolver::solve ( Problem *const &  problem_pt,
DoubleVector result 
)
pure virtual

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.

Implemented in oomph::TrilinosAztecOOSolver, oomph::HypreSolver, oomph::GMRESBlockPreconditioner, oomph::GMRESBlockPreconditioner, oomph::HelmholtzFGMRESMG< MATRIX >, oomph::HelmholtzGMRESMG< MATRIX >, oomph::ComplexGMRES< MATRIX >, oomph::ComplexDampedJacobi< MATRIX >, oomph::MumpsSolver, oomph::SuperLUSolver, oomph::FD_LU, oomph::DenseLU, oomph::AugmentedProblemGMRES, oomph::GMRES< MATRIX >, oomph::DampedJacobi< MATRIX >, oomph::GS< CRDoubleMatrix >, oomph::GS< MATRIX >, oomph::BiCGStab< MATRIX >, oomph::CG< MATRIX >, oomph::MGSolver< DIM >, oomph::HSL_MA42, oomph::BlockHopfLinearSolver, oomph::AugmentedBlockPitchForkLinearSolver, oomph::BlockPitchForkLinearSolver, and oomph::AugmentedBlockFoldLinearSolver.

Referenced by oomph::ProblemBasedShiftInvertOperator::apply(), oomph::Problem::calculate_continuation_derivatives(), oomph::FoldHandler::FoldHandler(), oomph::Problem::get_inverse_mass_matrix_times_residuals(), oomph::HopfHandler::HopfHandler(), oomph::Problem::newton_solve(), oomph::Problem::newton_solve_continuation(), run_it(), oomph::SolidICProblem::set_newmark_initial_condition_consistently(), oomph::HSL_MA42::solve(), oomph::BiCGStab< MATRIX >::solve(), oomph::GS< MATRIX >::solve(), oomph::GS< CRDoubleMatrix >::solve(), oomph::GMRES< MATRIX >::solve(), oomph::AugmentedProblemGMRES::solve(), oomph::FD_LU::solve(), oomph::ComplexGMRES< MATRIX >::solve(), oomph::HelmholtzGMRESMG< MATRIX >::solve(), oomph::DoubleMatrixBase::solve(), oomph::AugmentedBlockFoldLinearSolver::solve(), oomph::BlockPitchForkLinearSolver::solve(), oomph::AugmentedBlockPitchForkLinearSolver::solve(), oomph::BlockHopfLinearSolver::solve(), oomph::ARPACK::solve_eigenproblem(), and oomph::BlockHopfLinearSolver::solve_for_two_rhs().

◆ solve_transpose() [1/3]

virtual void oomph::LinearSolver::solve_transpose ( DoubleMatrixBase *const &  matrix_pt,
const DoubleVector rhs,
DoubleVector result 
)
inlinevirtual

Linear-algebra-type solver: Takes pointer to a matrix and rhs vector and returns the solution of the linear system.

Reimplemented in oomph::SuperLUSolver.

203  {
204  throw OomphLibError(
205  "DoubleVector based solve function not implemented for this solver",
208  }

References OOMPH_CURRENT_FUNCTION, and OOMPH_EXCEPTION_LOCATION.

◆ solve_transpose() [2/3]

virtual void oomph::LinearSolver::solve_transpose ( DoubleMatrixBase *const &  matrix_pt,
const Vector< double > &  rhs,
Vector< double > &  result 
)
inlinevirtual

Linear-algebra-type solver: Takes pointer to a matrix and rhs vector and returns the solution of the linear system.

215  {
216  throw OomphLibError(
217  "Vector<double> based solve function not implemented for this solver",
220  }

References OOMPH_CURRENT_FUNCTION, and OOMPH_EXCEPTION_LOCATION.

◆ solve_transpose() [3/3]

virtual void oomph::LinearSolver::solve_transpose ( Problem *const &  problem_pt,
DoubleVector result 
)
inlinevirtual

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).

Reimplemented in oomph::SuperLUSolver.

183  {
184  // Create an output stream
185  std::ostringstream error_message_stream;
186 
187  // Create the error message
188  error_message_stream << "The function to solve the transposed system has "
189  << "not yet been\nimplemented for this solver."
190  << std::endl;
191 
192  // Throw the error message
193  throw OomphLibError(error_message_stream.str(),
196  } // End of solve_transpose

References OOMPH_CURRENT_FUNCTION, and OOMPH_EXCEPTION_LOCATION.

Member Data Documentation

◆ Compute_gradient

bool oomph::LinearSolver::Compute_gradient
protected

flag that indicates whether the gradient required for the globally convergent Newton method should be computed or not

Referenced by disable_computation_of_gradient(), oomph::GMRES< MATRIX >::enable_computation_of_gradient(), oomph::SuperLUSolver::enable_computation_of_gradient(), oomph::SuperLUSolver::solve(), and oomph::SuperLUSolver::solve_transpose().

◆ Doc_time

◆ Enable_resolve

◆ Gradient_for_glob_conv_newton_solve

DoubleVector oomph::LinearSolver::Gradient_for_glob_conv_newton_solve
protected

DoubleVector storing the gradient for the globally convergent Newton method

Referenced by get_gradient(), reset_gradient(), oomph::SuperLUSolver::solve(), and oomph::SuperLUSolver::solve_transpose().

◆ Gradient_has_been_computed

bool oomph::LinearSolver::Gradient_has_been_computed
protected

flag that indicates whether the gradient was computed or not

Referenced by get_gradient(), oomph::SuperLUSolver::solve(), and oomph::SuperLUSolver::solve_transpose().


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