|
IBAMR
IBAMR version 0.19.
|
Class PETScNewtonKrylovSolver provides a NewtonKrylovSolver interface for a PETSc inexact Newton-Krylov iterative nonlinear solver (SNES). More...
#include <ibtk/PETScNewtonKrylovSolver.h>

Public Member Functions | |
| PETScNewtonKrylovSolver (std::string object_name, SAMRAI::tbox::Pointer< SAMRAI::tbox::Database > input_db, std::string default_options_prefix, MPI_Comm petsc_comm=PETSC_COMM_WORLD) | |
| Constructor for a concrete NewtonKrylovSolver that employs the PETSc SNES solver framework. More... | |
| PETScNewtonKrylovSolver (std::string object_name, SNES petsc_snes) | |
| Constructor for a concrete NewtonKrylovSolver that acts as a "wrapper" for a supplied PETSc SNES object. More... | |
| ~PETScNewtonKrylovSolver () | |
| Destructor. More... | |
| void | setOptionsPrefix (const std::string &options_prefix) |
| Set the options prefix used by this PETSc solver object. More... | |
| void | setHierarchyMathOps (SAMRAI::tbox::Pointer< HierarchyMathOps > hier_math_ops) override |
| Set the HierarchyMathOps object used by the solver. More... | |
General-purpose solver functionality. | |
| void | setHomogeneousBc (bool homogeneous_bc) override |
| Set whether the solver should use homogeneous boundary conditions. More... | |
| void | setSolutionTime (double solution_time) override |
| Set the time at which the solution is to be evaluated. More... | |
| void | setTimeInterval (double current_time, double new_time) override |
| Set the current time interval. More... | |
Functions to access solver parameters. | |
| virtual void | setMaxEvaluations (int max_evaluations) |
| Set the maximum number of function evaluations to use per solve. More... | |
| virtual int | getMaxEvaluations () const |
| Get the maximum number of function evaluations to use per solve. More... | |
| virtual void | setSolutionTolerance (double solution_tol) |
| Set the tolerance in terms of the norm of the change in the solution between steps. More... | |
| virtual double | getSolutionTolerance () const |
| Get the tolerance in terms of the norm of the change in the solution between steps. More... | |
Functions to access the underlying PETSc objects. | |
| const SNES & | getPETScSNES () const |
| Get the PETSc SNES object. More... | |
General-purpose solver functionality. | |
| void | setHomogeneousBc (bool homogeneous_bc) override |
| Set whether the solver should use homogeneous boundary conditions. More... | |
| void | setSolutionTime (double solution_time) override |
| Set the time at which the solution is to be evaluated. More... | |
| void | setTimeInterval (double current_time, double new_time) override |
| Set the current time interval. More... | |
Newton-Krylov solver functionality. | |
| virtual SAMRAI::tbox::Pointer< GeneralOperator > | getOperator () const |
| Retrieve the nonlinear operator \(F[x]\) used by the solver. More... | |
| virtual SAMRAI::tbox::Pointer< JacobianOperator > | getJacobian () const |
| Retrieve the Jacobian operator \(J[x] = F'[x]\) used by the solver. More... | |
| virtual SAMRAI::tbox::Pointer< KrylovLinearSolver > | getLinearSolver () const |
| Retrieve the Krylov linear solver used in computing Newton step directions. More... | |
Functions to access solver parameters. | |
| virtual void | setMaxEvaluations (int max_evaluations) |
| Set the maximum number of function evaluations to use per solve. More... | |
| virtual int | getMaxEvaluations () const |
| Get the maximum number of function evaluations to use per solve. More... | |
| virtual void | setSolutionTolerance (double solution_tol) |
| Set the tolerance in terms of the norm of the change in the solution between steps. More... | |
| virtual double | getSolutionTolerance () const |
| Get the tolerance in terms of the norm of the change in the solution between steps. More... | |
Static Public Member Functions | |
| static SAMRAI::tbox::Pointer< NewtonKrylovSolver > | allocate_solver (const std::string &object_name, SAMRAI::tbox::Pointer< SAMRAI::tbox::Database > input_db, const std::string &default_options_prefix) |
| Static function to construct a PETScNewtonKrylovSolver. More... | |
Functions to access solver parameters. | |
| virtual void | setMaxIterations (int max_iterations) |
| Set the maximum number of nonlinear iterations to use per solve. More... | |
| virtual int | getMaxIterations () const |
| Get the maximum number of nonlinear iterations to use per solve. More... | |
| virtual void | setAbsoluteTolerance (double abs_residual_tol) |
| Set the absolute residual tolerance for convergence. More... | |
| virtual double | getAbsoluteTolerance () const |
| Get the absolute residual tolerance for convergence. More... | |
| virtual void | setRelativeTolerance (double rel_residual_tol) |
| Set the relative residual tolerance for convergence. More... | |
| virtual double | getRelativeTolerance () const |
| Get the relative residual tolerance for convergence. More... | |
Logging functions. | |
| virtual void | setLoggingEnabled (bool enable_logging=true) |
| Enable or disable logging. More... | |
| virtual bool | getLoggingEnabled () const |
| Determine whether logging is enabled or disabled. More... | |
| virtual void | printClassData (std::ostream &stream) |
| Print class data to stream. More... | |
| 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 |
| void | init (const std::string &object_name, bool homogeneous_bc) |
| virtual void | initSpecialized (const std::string &object_name, bool homogeneous_bc) |
Functions to access data on the most recent solve. | |
| virtual int | getNumLinearIterations () const |
| Return the number of linear iterations from the most recent nonlinear solve. More... | |
| 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 |
Newton-Krylov solver functionality. | |
| void | setOperator (SAMRAI::tbox::Pointer< GeneralOperator > op) override |
| Set the nonlinear operator \(F[x]\) used by the solver. More... | |
| SAMRAI::tbox::Pointer< SAMRAI::solv::SAMRAIVectorReal< NDIM, double > > | getSolutionVector () const override |
| Return the vector in which the approximate solution is stored. More... | |
| SAMRAI::tbox::Pointer< SAMRAI::solv::SAMRAIVectorReal< NDIM, double > > | getFunctionVector () const override |
| Return the vector in which the nonlinear function evaluation is stored. More... | |
| void | setJacobian (SAMRAI::tbox::Pointer< JacobianOperator > J) override |
| Set the Jacobian operator \(J[x] = F'[x]\) used by the solver. More... | |
| bool | solveSystem (SAMRAI::solv::SAMRAIVectorReal< NDIM, double > &x, SAMRAI::solv::SAMRAIVectorReal< NDIM, double > &b) override |
| Solve the system \(F[x]=b\) for \(x\). More... | |
| void | initializeSolverState (const SAMRAI::solv::SAMRAIVectorReal< NDIM, double > &x, const SAMRAI::solv::SAMRAIVectorReal< NDIM, double > &b) override |
| Compute hierarchy dependent data required for solving \(F[x]=b\). More... | |
| void | deallocateSolverState () override |
| Remove all hierarchy dependent data allocated by initializeSolverState(). More... | |
| PETScNewtonKrylovSolver ()=delete | |
| Default constructor. More... | |
| PETScNewtonKrylovSolver (const PETScNewtonKrylovSolver &from)=delete | |
| Copy constructor. More... | |
| PETScNewtonKrylovSolver & | operator= (const PETScNewtonKrylovSolver &that)=delete |
| Assignment operator. More... | |
| void | common_ctor () |
| Common routine used by all class constructors. More... | |
| void | resetWrappedSNES (SNES &petsc_snes) |
| Reset the SNES wrapped by this solver class. More... | |
| void | resetSNESOptions () |
| Reset the values of the convergence tolerances for the PETSc SNES object. More... | |
| void | resetSNESFunction () |
| Reset the function for the PETSc SNES object. More... | |
| void | resetSNESJacobian () |
| Reset the Jacobian for the PETSc SNES object. More... | |
Static functions for use by PETSc SNES and MatShell objects. | |
| bool | d_reinitializing_solver = false |
| Vec | d_petsc_x = nullptr |
| Vec | d_petsc_b = nullptr |
| Vec | d_petsc_r = nullptr |
| std::string | d_options_prefix |
| MPI_Comm | d_petsc_comm = PETSC_COMM_WORLD |
| SNES | d_petsc_snes = nullptr |
| Mat | d_petsc_jac = nullptr |
| bool | d_managing_petsc_snes = true |
| bool | d_user_provided_function = false |
| bool | d_user_provided_jacobian = false |
| static PetscErrorCode | FormFunction_SAMRAI (SNES snes, Vec x, Vec f, void *p_ctx) |
| Evaluate f = F[x]. More... | |
| static PetscErrorCode | FormJacobian_SAMRAI (SNES snes, Vec x, Mat A, Mat B, void *p_ctx) |
| Setup F'[x]. More... | |
| static PetscErrorCode | MatVecMult_SAMRAI (Mat A, Vec x, Vec y) |
| Compute the matrix vector product y = Ax. More... | |
This solver class provides access to inexact Newton-Krylov methods, including line search and trust region Newton methods, using the PETSc SNES nonlinear solver interface. Note that solver configuration is typically done at runtime via command line options.
Sample parameters for initialization from database (and their default values):
options_prefix = "" // see setOptionsPrefix() rel_residual_tol = 1.0e-8 // see setRelativeTolerance() abs_residual_tol = 1.0e-50 // see setAbsoluteTolerance() solution_tol = 1.0e-8 // see setSolutionTolerance() max_iterations = 50 // see setMaxIterations() enable_logging = FALSE // see setLoggingEnabled()
PETSc is developed in the Mathematics and Computer Science (MCS) Division at Argonne National Laboratory (ANL). For more information about PETSc, see http://www.mcs.anl.gov/petsc.
| IBTK::PETScNewtonKrylovSolver::PETScNewtonKrylovSolver | ( | std::string | object_name, |
| SAMRAI::tbox::Pointer< SAMRAI::tbox::Database > | input_db, | ||
| std::string | default_options_prefix, | ||
| MPI_Comm | petsc_comm = PETSC_COMM_WORLD |
||
| ) |
| IBTK::PETScNewtonKrylovSolver::PETScNewtonKrylovSolver | ( | std::string | object_name, |
| SNES | petsc_snes | ||
| ) |
| object_name | Name of the solver |
| petsc_snes | PETSc SNES object |
| IBTK::PETScNewtonKrylovSolver::~PETScNewtonKrylovSolver | ( | ) |
|
privatedelete |
|
privatedelete |
| from | The value to copy to this object. |
|
inlinestatic |
| void IBTK::PETScNewtonKrylovSolver::setOptionsPrefix | ( | const std::string & | options_prefix | ) |
| const SNES& IBTK::PETScNewtonKrylovSolver::getPETScSNES | ( | ) | const |
|
overridevirtual |
Reimplemented from IBTK::NewtonKrylovSolver.
|
overridevirtual |
Implements IBTK::NewtonKrylovSolver.
|
overridevirtual |
Implements IBTK::NewtonKrylovSolver.
|
overridevirtual |
Reimplemented from IBTK::NewtonKrylovSolver.
|
overridevirtual |
Before calling solveSystem(), the form of the solution x and right-hand-side b vectors must be set properly by the user on all patch interiors on the specified range of levels in the patch hierarchy. The user is responsible for all data management for the quantities associated with the solution and right-hand-side vectors. In particular, patch data in these vectors must be allocated prior to calling this method.
| x | solution vector |
| b | right-hand-side vector |
Conditions on Parameters:
true if the solver converged to the specified tolerances, false otherwise Implements IBTK::GeneralSolver.
|
overridevirtual |
By default, the solveSystem() method computes some required hierarchy dependent data before solving and removes that data after the solve. For multiple solves that use the same hierarchy configuration, it is more efficient to:
Note that it is generally necessary to reinitialize the solver state when the hierarchy configuration changes.
When any operator objects have been registered with this class via setOperator() or setJacobian(), they are also initialized by this member function.
| x | solution vector |
| b | right-hand-side vector |
Conditions on Parameters:
Reimplemented from IBTK::GeneralSolver.
|
overridevirtual |
When any operator objects have been registered with this class via setOperator() or setJacobian(), they are also deallocated by this member function.
Reimplemented from IBTK::GeneralSolver.
|
privatedelete |
| that | The value to assign to this object. |
|
private |
|
private |
|
private |
|
private |
|
private |
|
staticprivate |
|
staticprivate |
|
staticprivate |
|
overridevirtualinherited |
Reimplemented from IBTK::GeneralSolver.
|
overridevirtualinherited |
Reimplemented from IBTK::GeneralSolver.
|
overridevirtualinherited |
Reimplemented from IBTK::GeneralSolver.
|
overridevirtualinherited |
Reimplemented from IBTK::GeneralSolver.
|
virtualinherited |
|
virtualinherited |
|
virtualinherited |
|
virtualinherited |
|
virtualinherited |
|
virtualinherited |
|
virtualinherited |
|
virtualinherited |
|
inherited |
|
virtualinherited |
|
virtualinherited |
|
virtualinherited |
|
virtualinherited |
|
virtualinherited |
|
virtualinherited |
|
virtualinherited |
Reimplemented in IBTK::PETScPCLSWrapper.
|
virtualinherited |
Reimplemented in IBTK::PETScPCLSWrapper.
|
virtualinherited |
Reimplemented in IBTK::PETScPCLSWrapper.
|
virtualinherited |
Reimplemented in IBTK::PETScPCLSWrapper.
|
virtualinherited |
Reimplemented in IBTK::PETScPCLSWrapper.
|
virtualinherited |
Reimplemented in IBTK::BGaussSeidelPreconditioner, IBTK::PETScPCLSWrapper, and IBTK::BJacobiPreconditioner.
|
virtualinherited |
Reimplemented in IBTK::BGaussSeidelPreconditioner, IBTK::PETScPCLSWrapper, and IBTK::BJacobiPreconditioner.
|
virtualinherited |
|
virtualinherited |
|
virtualinherited |
Reimplemented in IBTK::LinearSolver.
|
protectedinherited |
|
protectedvirtualinherited |
Reimplemented in IBTK::PoissonSolver.
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
1.8.17