|
IBAMR
An adaptive and distributed-memory parallel implementation of the immersed boundary (IB) method
|
Class NewtonKrylovSolver provides an abstract interface for the implementation of inexact Newton-Krylov solvers for nonlinear problems of the form ![$ F[x]=b $](form_227.png)
#include </home/runner/work/IBAMR/IBAMR/ibtk/include/ibtk/NewtonKrylovSolver.h>

Public Member Functions | |
| NewtonKrylovSolver () | |
| Default constructor. | |
| virtual | ~NewtonKrylovSolver ()=default |
| Empty virtual destructor. | |
| void | setHierarchyMathOps (SAMRAI::tbox::Pointer< HierarchyMathOps > hier_math_ops) override |
| Set the HierarchyMathOps object used by the solver. | |
General-purpose solver functionality. | |
| void | setHomogeneousBc (bool homogeneous_bc) override |
| Set whether the solver should use homogeneous boundary conditions. | |
| void | setSolutionTime (double solution_time) override |
| Set the time at which the solution is to be evaluated. | |
| void | setTimeInterval (double current_time, double new_time) override |
| Set the current time interval. | |
Newton-Krylov solver functionality. | |
| virtual void | setOperator (SAMRAI::tbox::Pointer< GeneralOperator > op) |
Set the nonlinear operator ![]() | |
| virtual SAMRAI::tbox::Pointer< GeneralOperator > | getOperator () const |
Retrieve the nonlinear operator ![]() | |
| virtual SAMRAI::tbox::Pointer< SAMRAI::solv::SAMRAIVectorReal< NDIM, double > > | getSolutionVector () const =0 |
| Return the vector in which the approximate solution is stored. | |
| virtual SAMRAI::tbox::Pointer< SAMRAI::solv::SAMRAIVectorReal< NDIM, double > > | getFunctionVector () const =0 |
| Return the vector in which the nonlinear function evaluation is stored. | |
| virtual void | setJacobian (SAMRAI::tbox::Pointer< JacobianOperator > J) |
Set the Jacobian operator ![]() | |
| virtual SAMRAI::tbox::Pointer< JacobianOperator > | getJacobian () const |
Retrieve the Jacobian operator ![]() | |
| virtual SAMRAI::tbox::Pointer< KrylovLinearSolver > | getLinearSolver () const |
| Retrieve the Krylov linear solver used in computing Newton step directions. | |
Functions to access solver parameters. | |
| virtual void | setMaxEvaluations (int max_evaluations) |
| Set the maximum number of function evaluations to use per solve. | |
| virtual int | getMaxEvaluations () const |
| Get the maximum number of function evaluations to use per solve. | |
| virtual void | setSolutionTolerance (double solution_tol) |
| Set the tolerance in terms of the norm of the change in the solution between steps. | |
| virtual double | getSolutionTolerance () const |
| Get the tolerance in terms of the norm of the change in the solution between steps. | |
Public Member Functions inherited from IBTK::GeneralSolver | |
| GeneralSolver ()=default | |
| Constructor. | |
| virtual | ~GeneralSolver ()=default |
| Empty virtual destructor. | |
| const std::string & | getName () const |
| Return the object name. | |
| virtual bool | getIsInitialized () const |
| Return whether the operator is initialized. | |
| virtual bool | getHomogeneousBc () const |
| Return whether the solver is using homogeneous boundary conditions. | |
| virtual double | getSolutionTime () const |
| Get the time at which the solution is being evaluated. | |
| virtual std::pair< double, double > | getTimeInterval () const |
| Get the current time interval. | |
| virtual double | getDt () const |
| Get the current time step size. | |
| virtual SAMRAI::tbox::Pointer< HierarchyMathOps > | getHierarchyMathOps () const |
| Get the HierarchyMathOps object used by the solver. | |
| virtual bool | solveSystem (SAMRAI::solv::SAMRAIVectorReal< NDIM, double > &x, SAMRAI::solv::SAMRAIVectorReal< NDIM, double > &b)=0 |
| Solve the system of equations. | |
| virtual void | initializeSolverState (const SAMRAI::solv::SAMRAIVectorReal< NDIM, double > &x, const SAMRAI::solv::SAMRAIVectorReal< NDIM, double > &b) |
Compute hierarchy dependent data required for solving ![]() | |
| virtual void | deallocateSolverState () |
| Remove all hierarchy dependent data allocated by initializeSolverState(). | |
| virtual void | setMaxIterations (int max_iterations) |
| Set the maximum number of nonlinear iterations to use per solve. | |
| virtual int | getMaxIterations () const |
| Get the maximum number of nonlinear iterations to use per solve. | |
| virtual void | setAbsoluteTolerance (double abs_residual_tol) |
| Set the absolute residual tolerance for convergence. | |
| virtual double | getAbsoluteTolerance () const |
| Get the absolute residual tolerance for convergence. | |
| virtual void | setRelativeTolerance (double rel_residual_tol) |
| Set the relative residual tolerance for convergence. | |
| virtual double | getRelativeTolerance () const |
| Get the relative residual tolerance for convergence. | |
| virtual int | getNumIterations () const |
| Return the iteration count from the most recent solve. | |
| virtual double | getResidualNorm () const |
| Return the residual norm from the most recent iteration. | |
| virtual void | setLoggingEnabled (bool enable_logging=true) |
| Enable or disable logging. | |
| virtual bool | getLoggingEnabled () const |
| Determine whether logging is enabled or disabled. | |
| virtual void | printClassData (std::ostream &stream) |
| Print class data to stream. | |
Functions to access data on the most recent solve. | |
| SAMRAI::tbox::Pointer< GeneralOperator > | d_F |
| SAMRAI::tbox::Pointer< JacobianOperator > | d_J |
| SAMRAI::tbox::Pointer< KrylovLinearSolver > | d_krylov_solver |
| SAMRAI::tbox::Pointer< SAMRAI::solv::SAMRAIVectorReal< NDIM, double > > | d_x |
| SAMRAI::tbox::Pointer< SAMRAI::solv::SAMRAIVectorReal< NDIM, double > > | d_b |
| SAMRAI::tbox::Pointer< SAMRAI::solv::SAMRAIVectorReal< NDIM, double > > | d_r |
| int | d_max_evaluations = 10000 |
| double | d_solution_tol = 1.0e-8 |
| int | d_current_linear_iterations = 0 |
| virtual int | getNumLinearIterations () const |
| Return the number of linear iterations from the most recent nonlinear solve. | |
Additional Inherited Members | |
Protected Member Functions inherited from IBTK::GeneralSolver | |
| void | init (const std::string &object_name, bool homogeneous_bc) |
| virtual void | initSpecialized (const std::string &object_name, bool homogeneous_bc) |
Protected Attributes inherited from IBTK::GeneralSolver | |
| std::string | d_object_name = "unitialized" |
| bool | d_is_initialized = false |
| bool | d_homogeneous_bc = false |
| double | d_solution_time = std::numeric_limits<double>::quiet_NaN() |
| double | d_current_time = std::numeric_limits<double>::quiet_NaN() |
| double | d_new_time = std::numeric_limits<double>::quiet_NaN() |
| double | d_rel_residual_tol = 0.0 |
| double | d_abs_residual_tol = 0.0 |
| int | d_max_iterations = 100 |
| int | d_current_iterations = 0 |
| double | d_current_residual_norm = std::numeric_limits<double>::quiet_NaN() |
| SAMRAI::tbox::Pointer< HierarchyMathOps > | d_hier_math_ops |
| bool | d_hier_math_ops_external = false |
| bool | d_enable_logging = false |
Class NewtonKrylovSolver provides an abstract interface for the implementation of inexact Newton-Krylov solvers for nonlinear problems of the form ![$ F[x]=b $](form_227.png)
|
pure virtual |
Return the vector in which the nonlinear function evaluation is stored.
Implemented in IBTK::PETScNewtonKrylovSolver.
|
pure virtual |
Return the vector in which the approximate solution is stored.
Implemented in IBTK::PETScNewtonKrylovSolver.
|
overridevirtual |
Set the HierarchyMathOps object used by the solver.
Reimplemented from IBTK::GeneralSolver.
Set whether the solver should use homogeneous boundary conditions.
Reimplemented from IBTK::GeneralSolver.
|
virtual |
Set the Jacobian operator ![$J[x] = F'[x]$](form_229.png)
Reimplemented in IBTK::PETScNewtonKrylovSolver.
|
virtual |
Set the nonlinear operator ![$F[x]$](form_228.png)
Reimplemented in IBTK::PETScNewtonKrylovSolver.
Set the time at which the solution is to be evaluated.
Reimplemented from IBTK::GeneralSolver.
|
overridevirtual |
Set the current time interval.
Reimplemented from IBTK::GeneralSolver.