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). More...
#include <linear_solver.h>
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 () |
Enable resolve (i.e. store matrix and/or LU decomposition, say) Virtual so it can be overloaded to perform additional tasks. More... | |
virtual void | disable_resolve () |
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. More... | |
virtual void | solve (Problem *const &problem_pt, DoubleVector &result)=0 |
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. More... | |
virtual void | solve (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. More... | |
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. More... | |
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). More... | |
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. More... | |
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. More... | |
virtual void | resolve (const DoubleVector &rhs, DoubleVector &result) |
Resolve the system defined by the last assembled jacobian and the rhs vector. Solution is returned in the vector result. (broken virtual) More... | |
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) More... | |
virtual void | clean_up_memory () |
Empty virtual function that can be overloaded in specific linear solvers to clean up any memory that may have been allocated (e.g. when preparing for a re-solve). More... | |
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) More... | |
virtual double | linear_solver_solution_time () const |
return the time taken to solve the linear system (needs to be overloaded for each linear solver) More... | |
virtual void | enable_computation_of_gradient () |
function to enable the computation of the gradient required for the globally convergent Newton method More... | |
void | disable_computation_of_gradient () |
function to disable the computation of the gradient required for the globally convergent Newton method More... | |
void | reset_gradient () |
function to reset the size of the gradient before each Newton solve More... | |
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... | |
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... | |
Protected Attributes | |
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. More... | |
bool | Doc_time |
Boolean flag that indicates whether the time taken. More... | |
bool | Compute_gradient |
flag that indicates whether the gradient required for the globally convergent Newton method should be computed or not More... | |
bool | Gradient_has_been_computed |
flag that indicates whether the gradient was computed or not More... | |
DoubleVector | Gradient_for_glob_conv_newton_solve |
DoubleVector storing the gradient for the globally convergent Newton method. More... | |
Additional Inherited Members | |
Protected Member Functions inherited from oomph::DistributableLinearAlgebraObject | |
void | clear_distribution () |
clear the distribution of this distributable linear algebra object More... | |
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).
Definition at line 67 of file linear_solver.h.
|
inline |
Empty constructor, initialise the member data.
Definition at line 92 of file linear_solver.h.
|
delete |
Broken copy constructor.
|
inlinevirtual |
Empty virtual destructor.
Definition at line 107 of file linear_solver.h.
|
inlinevirtual |
Empty virtual function that can be overloaded in specific linear solvers to clean up any memory that may have been allocated (e.g. when preparing for a re-solve).
Reimplemented in oomph::GMRESBlockPreconditioner, oomph::GMRESBlockPreconditioner, oomph::HelmholtzGMRESMG< MATRIX >, oomph::ComplexGMRES< MATRIX >, oomph::ComplexDampedJacobi< MATRIX >, oomph::TrilinosAztecOOSolver, oomph::MumpsSolver, oomph::SuperLUSolver, oomph::DenseLU, oomph::AugmentedProblemGMRES, oomph::GMRES< MATRIX >, oomph::DampedJacobi< MATRIX >, oomph::GS< CRDoubleMatrix >, oomph::GS< MATRIX >, oomph::BiCGStab< MATRIX >, oomph::CG< MATRIX >, oomph::HypreSolver, oomph::MGPreconditioner< DIM >, oomph::MGSolver< DIM >, and oomph::HSL_MA42.
Definition at line 256 of file linear_solver.h.
Referenced by oomph::CRDoubleMatrix::clear().
|
inline |
function to disable the computation of the gradient required for the globally convergent Newton method
Definition at line 293 of file linear_solver.h.
References Compute_gradient.
|
inline |
Disable documentation of solve times.
Definition at line 116 of file linear_solver.h.
References Doc_time.
Referenced by oomph::AdjointProblemBasedShiftInvertOperator::AdjointProblemBasedShiftInvertOperator(), oomph::NewMumpsPreconditioner::disable_doc_time(), oomph::HelmholtzMGPreconditioner< DIM >::disable_smoother_and_superlu_doc_time(), oomph::ProblemBasedShiftInvertOperator::ProblemBasedShiftInvertOperator(), oomph::ProjectionProblem< PROJECTABLE_ELEMENT >::project(), oomph::ARPACK::solve_eigenproblem_legacy(), and oomph::SuperLUPreconditioner::SuperLUPreconditioner().
|
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.
Definition at line 144 of file linear_solver.h.
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(), oomph::AugmentedBlockFoldLinearSolver::solve(), oomph::BlockPitchForkLinearSolver::solve(), oomph::AugmentedBlockPitchForkLinearSolver::solve(), oomph::BlockHopfLinearSolver::solve(), oomph::ARPACK::solve_eigenproblem_legacy(), and oomph::BlockHopfLinearSolver::solve_for_two_rhs().
|
inlinevirtual |
function to enable the computation of the gradient required for the globally convergent Newton method
Reimplemented in oomph::SuperLUSolver, and oomph::GMRES< MATRIX >.
Definition at line 282 of file linear_solver.h.
Referenced by oomph::Problem::newton_solve().
|
inline |
Enable documentation of solve times.
Definition at line 110 of file linear_solver.h.
References Doc_time.
Referenced by oomph::NewMumpsPreconditioner::enable_doc_time(), and oomph::ProjectionProblem< PROJECTABLE_ELEMENT >::project().
|
inlinevirtual |
Enable resolve (i.e. store matrix and/or LU decomposition, say) Virtual so it can be overloaded to perform additional tasks.
Definition at line 135 of file linear_solver.h.
References Enable_resolve.
Referenced by oomph::ProblemBasedShiftInvertOperator::apply(), oomph::AdjointProblemBasedShiftInvertOperator::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(), oomph::AugmentedBlockFoldLinearSolver::solve(), oomph::BlockPitchForkLinearSolver::solve(), oomph::AugmentedBlockPitchForkLinearSolver::solve(), oomph::BlockHopfLinearSolver::solve(), oomph::ARPACK::solve_eigenproblem_legacy(), and oomph::BlockHopfLinearSolver::solve_for_two_rhs().
|
inline |
function to access the gradient, provided it has been computed
Definition at line 306 of file linear_solver.h.
References Gradient_for_glob_conv_newton_solve, and Gradient_has_been_computed.
Referenced by oomph::Problem::newton_solve().
|
inline |
Is documentation of solve times enabled?
Definition at line 122 of file linear_solver.h.
References Doc_time.
Referenced by oomph::ProjectionProblem< PROJECTABLE_ELEMENT >::project().
|
inline |
Boolean flag indicating if resolves are enabled.
Definition at line 128 of file linear_solver.h.
References Enable_resolve.
Referenced by oomph::Problem::calculate_continuation_derivatives(), oomph::FoldHandler::FoldHandler(), oomph::HopfHandler::HopfHandler(), and oomph::Problem::newton_solve_continuation().
|
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.
Definition at line 260 of file linear_solver.h.
|
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.
Definition at line 271 of file linear_solver.h.
|
delete |
Broken assignment operator.
|
inline |
function to reset the size of the gradient before each Newton solve
Definition at line 300 of file linear_solver.h.
References oomph::DoubleVector::clear(), and Gradient_for_glob_conv_newton_solve.
Referenced by oomph::Problem::newton_solve().
|
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.
Definition at line 225 of file linear_solver.h.
Referenced by oomph::ProblemBasedShiftInvertOperator::apply(), oomph::AdjointProblemBasedShiftInvertOperator::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(), oomph::AugmentedBlockFoldLinearSolver::solve(), oomph::BlockPitchForkLinearSolver::solve(), oomph::AugmentedBlockPitchForkLinearSolver::solve(), oomph::BlockHopfLinearSolver::solve(), oomph::ARPACK::solve_eigenproblem_legacy(), and oomph::BlockHopfLinearSolver::solve_for_two_rhs().
|
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.
Definition at line 236 of file linear_solver.h.
|
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::TrilinosAztecOOSolver, oomph::AugmentedProblemGMRES, oomph::GMRES< MATRIX >, oomph::DampedJacobi< MATRIX >, oomph::GS< CRDoubleMatrix >, oomph::GS< MATRIX >, oomph::BiCGStab< MATRIX >, oomph::CG< MATRIX >, oomph::HypreSolver, oomph::MumpsSolver, oomph::SuperLUSolver, oomph::FD_LU, oomph::DenseLU, oomph::HSL_MA42, oomph::BlockHopfLinearSolver, oomph::AugmentedBlockPitchForkLinearSolver, oomph::BlockPitchForkLinearSolver, and oomph::AugmentedBlockFoldLinearSolver.
Definition at line 156 of file linear_solver.h.
|
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.
Definition at line 168 of file linear_solver.h.
|
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::AdjointProblemBasedShiftInvertOperator::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::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_legacy(), and oomph::BlockHopfLinearSolver::solve_for_two_rhs().
|
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.
Definition at line 200 of file linear_solver.h.
|
inlinevirtual |
Linear-algebra-type solver: Takes pointer to a matrix and rhs vector and returns the solution of the linear system.
Definition at line 212 of file linear_solver.h.
|
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.
Definition at line 181 of file linear_solver.h.
|
protected |
flag that indicates whether the gradient required for the globally convergent Newton method should be computed or not
Definition at line 81 of file linear_solver.h.
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().
|
protected |
Boolean flag that indicates whether the time taken.
Definition at line 77 of file linear_solver.h.
Referenced by oomph::MumpsSolver::backsub(), oomph::DenseLU::DenseLU(), disable_doc_time(), oomph::MGSolver< DIM >::disable_output(), oomph::MGSolver< DIM >::disable_v_cycle_output(), enable_doc_time(), oomph::MGSolver< DIM >::enable_output(), oomph::MGSolver< DIM >::enable_v_cycle_output(), oomph::MumpsSolver::factorise(), oomph::HypreSolver::HypreSolver(), is_doc_time_enabled(), oomph::GMRESBlockPreconditioner::preconditioner_solve(), oomph::SuperLUSolver::resolve(), oomph::MumpsSolver::resolve(), oomph::HypreSolver::resolve(), oomph::TrilinosAztecOOSolver::resolve(), oomph::SuperLUSolver::resolve_transpose(), oomph::MumpsSolver::solve(), oomph::GS< CRDoubleMatrix >::solve(), oomph::AugmentedProblemGMRES::solve(), oomph::DenseLU::solve(), oomph::FD_LU::solve(), oomph::SuperLUSolver::solve(), oomph::HelmholtzGMRESMG< MATRIX >::solve(), oomph::HelmholtzFGMRESMG< MATRIX >::solve(), oomph::HypreSolver::solve(), oomph::TrilinosAztecOOSolver::solve(), oomph::GS< CRDoubleMatrix >::solve_helper(), oomph::AugmentedProblemGMRES::solve_helper(), oomph::SuperLUSolver::solve_transpose(), oomph::TrilinosAztecOOSolver::solve_using_AztecOO(), and oomph::TrilinosAztecOOSolver::solver_setup().
|
protected |
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.
Definition at line 73 of file linear_solver.h.
Referenced by oomph::HypreSolver::disable_resolve(), disable_resolve(), oomph::TrilinosAztecOOSolver::disable_resolve(), enable_resolve(), is_resolve_enabled(), oomph::DenseLU::solve(), oomph::SuperLUSolver::solve(), oomph::MumpsSolver::solve(), oomph::CG< MATRIX >::solve(), oomph::BiCGStab< MATRIX >::solve(), oomph::GS< MATRIX >::solve(), oomph::DampedJacobi< MATRIX >::solve(), oomph::GMRES< MATRIX >::solve(), oomph::TrilinosAztecOOSolver::solve(), oomph::AugmentedBlockFoldLinearSolver::solve(), oomph::BlockPitchForkLinearSolver::solve(), oomph::AugmentedBlockPitchForkLinearSolver::solve(), oomph::BlockHopfLinearSolver::solve(), oomph::HSL_MA42::solve(), oomph::GS< CRDoubleMatrix >::solve(), oomph::AugmentedProblemGMRES::solve(), oomph::HelmholtzGMRESMG< MATRIX >::solve(), oomph::HelmholtzFGMRESMG< MATRIX >::solve(), oomph::HypreSolver::solve(), oomph::HSL_MA42::solve_for_one_dof(), oomph::BlockHopfLinearSolver::solve_for_two_rhs(), and oomph::SuperLUSolver::solve_transpose().
|
protected |
DoubleVector storing the gradient for the globally convergent Newton method.
Definition at line 88 of file linear_solver.h.
Referenced by get_gradient(), reset_gradient(), oomph::SuperLUSolver::solve(), and oomph::SuperLUSolver::solve_transpose().
|
protected |
flag that indicates whether the gradient was computed or not
Definition at line 84 of file linear_solver.h.
Referenced by get_gradient(), oomph::SuperLUSolver::solve(), and oomph::SuperLUSolver::solve_transpose().