IBAMR
An adaptive and distributed-memory parallel implementation of the immersed boundary (IB) method
Public Member Functions | Static Public Member Functions | List of all members
IBTK::PETScNewtonKrylovSolver Class Reference

Class PETScNewtonKrylovSolver provides a NewtonKrylovSolver interface for a PETSc inexact Newton-Krylov iterative nonlinear solver (SNES). More...

#include </home/runner/work/IBAMR/IBAMR/ibtk/include/ibtk/PETScNewtonKrylovSolver.h>

Inheritance diagram for IBTK::PETScNewtonKrylovSolver:
Inheritance graph
[legend]

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.
 
 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.
 
void setOptionsPrefix (const std::string &options_prefix)
 Set the options prefix used by this PETSc solver object.
 
Functions to access the underlying PETSc objects.
const SNES & getPETScSNES () const
 Get the PETSc SNES object.
 
- Public Member Functions inherited from IBTK::NewtonKrylovSolver
 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.
 
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.
 
virtual SAMRAI::tbox::Pointer< GeneralOperatorgetOperator () const
 Retrieve the nonlinear operator $F[x]$ used by the solver.
 
virtual SAMRAI::tbox::Pointer< JacobianOperatorgetJacobian () const
 Retrieve the Jacobian operator $J[x] = F'[x]$ used by the solver.
 
virtual SAMRAI::tbox::Pointer< KrylovLinearSolvergetLinearSolver () const
 Retrieve the Krylov linear solver used in computing Newton step directions.
 
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.
 
virtual int getNumLinearIterations () const
 Return the number of linear iterations from the most recent nonlinear solve.
 
- Public Member Functions inherited from IBTK::GeneralSolver
 GeneralSolver ()=default
 Constructor.
 
virtual ~GeneralSolver ()=default
 Empty virtual destructor.
 
const std::stringgetName () 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< HierarchyMathOpsgetHierarchyMathOps () const
 Get the HierarchyMathOps object used by the solver.
 
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.
 

Static Public Member Functions

static SAMRAI::tbox::Pointer< NewtonKrylovSolverallocate_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.
 

Newton-Krylov solver functionality.

void setOperator (SAMRAI::tbox::Pointer< GeneralOperator > op) override
 Set the nonlinear operator $F[x]$ used by the solver.
 
SAMRAI::tbox::Pointer< SAMRAI::solv::SAMRAIVectorReal< NDIM, double > > getSolutionVector () const override
 Return the vector in which the approximate solution is stored.
 
SAMRAI::tbox::Pointer< SAMRAI::solv::SAMRAIVectorReal< NDIM, double > > getFunctionVector () const override
 Return the vector in which the nonlinear function evaluation is stored.
 
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...
 

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::NewtonKrylovSolver
SAMRAI::tbox::Pointer< GeneralOperatord_F
 
SAMRAI::tbox::Pointer< JacobianOperatord_J
 
SAMRAI::tbox::Pointer< KrylovLinearSolverd_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
 
- 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< HierarchyMathOpsd_hier_math_ops
 
bool d_hier_math_ops_external = false
 
bool d_enable_logging = false
 

Detailed Description

Class PETScNewtonKrylovSolver provides a NewtonKrylovSolver interface for a PETSc inexact Newton-Krylov iterative nonlinear solver (SNES).

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.

Note
  • Users have the option of supplying the PETScNewtonKrylovSolver class constructor with a SNES object. It is important to emphasize that doing so is optional. When a SNES object is provided to the class constructor, the PETScNewtonKrylovSolver class acts as a wrapper for the supplied SNES object.
  • The functionality of the PETScNewtonKrylovSolver class is essentially identical regardless as to whether a PETSc SNES object is provided to the class constructor. The main exception is that when an external SNES object is provided to the class constructor, memory management of that object is NOT handled by the PETScNewtonKrylovSolver. In particular, it is the caller's responsibility to ensure that the supplied SNES object is properly destroyed via SNESDestroy().

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.

Constructor & Destructor Documentation

◆ PETScNewtonKrylovSolver()

IBTK::PETScNewtonKrylovSolver::PETScNewtonKrylovSolver ( std::string  object_name,
SNES  petsc_snes 
)

Constructor for a concrete NewtonKrylovSolver that acts as a "wrapper" for a supplied PETSc SNES object.

Parameters
object_nameName of the solver
petsc_snesPETSc SNES object
Note
This constructor initializes a PETScNewtonKrylovSolver object that acts as a "wrapper" for the provided SNES object. Note that memory management of the provided SNES object is NOT handled by this class.

Member Function Documentation

◆ deallocateSolverState()

void IBTK::PETScNewtonKrylovSolver::deallocateSolverState ( )
overridevirtual

Remove all hierarchy dependent data allocated by initializeSolverState().

When any operator objects have been registered with this class via setOperator() or setJacobian(), they are also deallocated by this member function.

Note
It is safe to call deallocateSolverState() when the solver state is already deallocated.
See also
initializeSolverState

Reimplemented from IBTK::GeneralSolver.

◆ initializeSolverState()

void IBTK::PETScNewtonKrylovSolver::initializeSolverState ( const SAMRAI::solv::SAMRAIVectorReal< NDIM, double > &  x,
const SAMRAI::solv::SAMRAIVectorReal< NDIM, double > &  b 
)
overridevirtual

Compute hierarchy dependent data required for solving $F[x]=b$.

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:

  1. initialize the hierarchy-dependent data required by the solver via initializeSolverState(),
  2. solve the system one or more times via solveSystem(), and
  3. remove the hierarchy-dependent data via deallocateSolverState().

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.

Parameters
xsolution vector
bright-hand-side vector

Conditions on Parameters:

  • vectors x and b must have same patch hierarchy
  • vectors x and b must have same structure, depth, etc.
Note
The vector arguments for solveSystem() need not match those for initializeSolverState(). However, there must be a certain degree of similarity, including:
  • hierarchy configuration (hierarchy pointer and range of levels)
  • number, type and alignment of vector component data
  • ghost cell widths of data in the solution x and right-hand-side b vectors
Note
It is safe to call initializeSolverState() when the state is already initialized. In this case, the solver state is first deallocated and then reinitialized.
See also
deallocateSolverState

Reimplemented from IBTK::GeneralSolver.

◆ setJacobian()

void IBTK::PETScNewtonKrylovSolver::setJacobian ( SAMRAI::tbox::Pointer< JacobianOperator J)
overridevirtual

Set the Jacobian operator $J[x] = F'[x]$ used by the solver.

Note
If a Jacobian object is not explicitly provided to the solver, a Jacobian-free inexact Newton-Krylov method is employed to approximate the action of the Jacobian.

Reimplemented from IBTK::NewtonKrylovSolver.

◆ solveSystem()

bool IBTK::PETScNewtonKrylovSolver::solveSystem ( SAMRAI::solv::SAMRAIVectorReal< NDIM, double > &  x,
SAMRAI::solv::SAMRAIVectorReal< NDIM, double > &  b 
)
overridevirtual

Solve the system $F[x]=b$ for $x$.

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.

Parameters
xsolution vector
bright-hand-side vector

Conditions on Parameters:

  • vectors x and b must have same patch hierarchy
  • vectors x and b must have same structure, depth, etc.
Note
The vector arguments for solveSystem() need not match those for initializeSolverState(). However, there must be a certain degree of similarity, including:
  • hierarchy configuration (hierarchy pointer and range of levels)
  • number, type and alignment of vector component data
  • ghost cell widths of data in the solution x and right-hand-side b vectors
Note
The solver need not be initialized prior to calling solveSystem(); however, see initializeSolverState() and deallocateSolverState() for opportunities to save overhead when performing multiple consecutive solves.
See also
initializeSolverState
deallocateSolverState
Returns
true if the solver converged to the specified tolerances, false otherwise

Implements IBTK::GeneralSolver.


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