IBAMR  IBAMR version 0.19.
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 <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. 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< GeneralOperatorgetOperator () const
 Retrieve the nonlinear operator \(F[x]\) used by the solver. More...
 
virtual SAMRAI::tbox::Pointer< JacobianOperatorgetJacobian () const
 Retrieve the Jacobian operator \(J[x] = F'[x]\) used by the solver. More...
 
virtual SAMRAI::tbox::Pointer< KrylovLinearSolvergetLinearSolver () 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< 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. 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< HierarchyMathOpsd_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< 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
 

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...
 
PETScNewtonKrylovSolveroperator= (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...
 

Detailed Description

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() [1/4]

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 
)

◆ PETScNewtonKrylovSolver() [2/4]

IBTK::PETScNewtonKrylovSolver::PETScNewtonKrylovSolver ( std::string  object_name,
SNES  petsc_snes 
)
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.

◆ ~PETScNewtonKrylovSolver()

IBTK::PETScNewtonKrylovSolver::~PETScNewtonKrylovSolver ( )

◆ PETScNewtonKrylovSolver() [3/4]

IBTK::PETScNewtonKrylovSolver::PETScNewtonKrylovSolver ( )
privatedelete
Note
This constructor is not implemented and should not be used.

◆ PETScNewtonKrylovSolver() [4/4]

IBTK::PETScNewtonKrylovSolver::PETScNewtonKrylovSolver ( const PETScNewtonKrylovSolver from)
privatedelete
Note
This constructor is not implemented and should not be used.
Parameters
fromThe value to copy to this object.

Member Function Documentation

◆ allocate_solver()

static SAMRAI::tbox::Pointer<NewtonKrylovSolver> IBTK::PETScNewtonKrylovSolver::allocate_solver ( const std::string &  object_name,
SAMRAI::tbox::Pointer< SAMRAI::tbox::Database input_db,
const std::string &  default_options_prefix 
)
inlinestatic

◆ setOptionsPrefix()

void IBTK::PETScNewtonKrylovSolver::setOptionsPrefix ( const std::string &  options_prefix)

◆ getPETScSNES()

const SNES& IBTK::PETScNewtonKrylovSolver::getPETScSNES ( ) const

◆ setOperator()

void IBTK::PETScNewtonKrylovSolver::setOperator ( SAMRAI::tbox::Pointer< GeneralOperator op)
overridevirtual

Reimplemented from IBTK::NewtonKrylovSolver.

◆ getSolutionVector()

SAMRAI::tbox::Pointer<SAMRAI::solv::SAMRAIVectorReal<NDIM, double> > IBTK::PETScNewtonKrylovSolver::getSolutionVector ( ) const
overridevirtual

◆ getFunctionVector()

SAMRAI::tbox::Pointer<SAMRAI::solv::SAMRAIVectorReal<NDIM, double> > IBTK::PETScNewtonKrylovSolver::getFunctionVector ( ) const
overridevirtual

◆ setJacobian()

void IBTK::PETScNewtonKrylovSolver::setJacobian ( SAMRAI::tbox::Pointer< JacobianOperator J)
overridevirtual
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

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.

◆ initializeSolverState()

void IBTK::PETScNewtonKrylovSolver::initializeSolverState ( const SAMRAI::solv::SAMRAIVectorReal< NDIM, double > &  x,
const SAMRAI::solv::SAMRAIVectorReal< NDIM, double > &  b 
)
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:

  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.

◆ deallocateSolverState()

void IBTK::PETScNewtonKrylovSolver::deallocateSolverState ( )
overridevirtual

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.

◆ operator=()

PETScNewtonKrylovSolver& IBTK::PETScNewtonKrylovSolver::operator= ( const PETScNewtonKrylovSolver that)
privatedelete
Note
This operator is not implemented and should not be used.
Parameters
thatThe value to assign to this object.
Returns
A reference to this object.

◆ common_ctor()

void IBTK::PETScNewtonKrylovSolver::common_ctor ( )
private

◆ resetWrappedSNES()

void IBTK::PETScNewtonKrylovSolver::resetWrappedSNES ( SNES &  petsc_snes)
private

◆ resetSNESOptions()

void IBTK::PETScNewtonKrylovSolver::resetSNESOptions ( )
private

◆ resetSNESFunction()

void IBTK::PETScNewtonKrylovSolver::resetSNESFunction ( )
private

◆ resetSNESJacobian()

void IBTK::PETScNewtonKrylovSolver::resetSNESJacobian ( )
private

◆ FormFunction_SAMRAI()

static PetscErrorCode IBTK::PETScNewtonKrylovSolver::FormFunction_SAMRAI ( SNES  snes,
Vec  x,
Vec  f,
void *  p_ctx 
)
staticprivate

◆ FormJacobian_SAMRAI()

static PetscErrorCode IBTK::PETScNewtonKrylovSolver::FormJacobian_SAMRAI ( SNES  snes,
Vec  x,
Mat  A,
Mat  B,
void *  p_ctx 
)
staticprivate

◆ MatVecMult_SAMRAI()

static PetscErrorCode IBTK::PETScNewtonKrylovSolver::MatVecMult_SAMRAI ( Mat  A,
Vec  x,
Vec  y 
)
staticprivate

◆ setHierarchyMathOps()

void IBTK::NewtonKrylovSolver::setHierarchyMathOps ( SAMRAI::tbox::Pointer< HierarchyMathOps hier_math_ops)
overridevirtualinherited

Reimplemented from IBTK::GeneralSolver.

◆ setHomogeneousBc()

void IBTK::NewtonKrylovSolver::setHomogeneousBc ( bool  homogeneous_bc)
overridevirtualinherited

Reimplemented from IBTK::GeneralSolver.

◆ setSolutionTime()

void IBTK::NewtonKrylovSolver::setSolutionTime ( double  solution_time)
overridevirtualinherited

Reimplemented from IBTK::GeneralSolver.

◆ setTimeInterval()

void IBTK::NewtonKrylovSolver::setTimeInterval ( double  current_time,
double  new_time 
)
overridevirtualinherited

Reimplemented from IBTK::GeneralSolver.

◆ getOperator()

virtual SAMRAI::tbox::Pointer<GeneralOperator> IBTK::NewtonKrylovSolver::getOperator ( ) const
virtualinherited

◆ getJacobian()

virtual SAMRAI::tbox::Pointer<JacobianOperator> IBTK::NewtonKrylovSolver::getJacobian ( ) const
virtualinherited

◆ getLinearSolver()

virtual SAMRAI::tbox::Pointer<KrylovLinearSolver> IBTK::NewtonKrylovSolver::getLinearSolver ( ) const
virtualinherited

◆ setMaxEvaluations()

virtual void IBTK::NewtonKrylovSolver::setMaxEvaluations ( int  max_evaluations)
virtualinherited

◆ getMaxEvaluations()

virtual int IBTK::NewtonKrylovSolver::getMaxEvaluations ( ) const
virtualinherited

◆ setSolutionTolerance()

virtual void IBTK::NewtonKrylovSolver::setSolutionTolerance ( double  solution_tol)
virtualinherited

◆ getSolutionTolerance()

virtual double IBTK::NewtonKrylovSolver::getSolutionTolerance ( ) const
virtualinherited

◆ getNumLinearIterations()

virtual int IBTK::NewtonKrylovSolver::getNumLinearIterations ( ) const
virtualinherited

◆ getName()

const std::string& IBTK::GeneralSolver::getName ( ) const
inherited

◆ getIsInitialized()

virtual bool IBTK::GeneralSolver::getIsInitialized ( ) const
virtualinherited

◆ getHomogeneousBc()

virtual bool IBTK::GeneralSolver::getHomogeneousBc ( ) const
virtualinherited

◆ getSolutionTime()

virtual double IBTK::GeneralSolver::getSolutionTime ( ) const
virtualinherited

◆ getTimeInterval()

virtual std::pair<double, double> IBTK::GeneralSolver::getTimeInterval ( ) const
virtualinherited

◆ getDt()

virtual double IBTK::GeneralSolver::getDt ( ) const
virtualinherited

◆ getHierarchyMathOps()

virtual SAMRAI::tbox::Pointer<HierarchyMathOps> IBTK::GeneralSolver::getHierarchyMathOps ( ) const
virtualinherited

◆ setMaxIterations()

virtual void IBTK::GeneralSolver::setMaxIterations ( int  max_iterations)
virtualinherited

◆ getMaxIterations()

virtual int IBTK::GeneralSolver::getMaxIterations ( ) const
virtualinherited

Reimplemented in IBTK::PETScPCLSWrapper.

◆ setAbsoluteTolerance()

virtual void IBTK::GeneralSolver::setAbsoluteTolerance ( double  abs_residual_tol)
virtualinherited

Reimplemented in IBTK::PETScPCLSWrapper.

◆ getAbsoluteTolerance()

virtual double IBTK::GeneralSolver::getAbsoluteTolerance ( ) const
virtualinherited

Reimplemented in IBTK::PETScPCLSWrapper.

◆ setRelativeTolerance()

virtual void IBTK::GeneralSolver::setRelativeTolerance ( double  rel_residual_tol)
virtualinherited

Reimplemented in IBTK::PETScPCLSWrapper.

◆ getRelativeTolerance()

virtual double IBTK::GeneralSolver::getRelativeTolerance ( ) const
virtualinherited

Reimplemented in IBTK::PETScPCLSWrapper.

◆ getNumIterations()

virtual int IBTK::GeneralSolver::getNumIterations ( ) const
virtualinherited

◆ getResidualNorm()

virtual double IBTK::GeneralSolver::getResidualNorm ( ) const
virtualinherited

◆ setLoggingEnabled()

virtual void IBTK::GeneralSolver::setLoggingEnabled ( bool  enable_logging = true)
virtualinherited

◆ getLoggingEnabled()

virtual bool IBTK::GeneralSolver::getLoggingEnabled ( ) const
virtualinherited

◆ printClassData()

virtual void IBTK::GeneralSolver::printClassData ( std::ostream &  stream)
virtualinherited

Reimplemented in IBTK::LinearSolver.

◆ init()

void IBTK::GeneralSolver::init ( const std::string &  object_name,
bool  homogeneous_bc 
)
protectedinherited

◆ initSpecialized()

virtual void IBTK::GeneralSolver::initSpecialized ( const std::string &  object_name,
bool  homogeneous_bc 
)
protectedvirtualinherited

Reimplemented in IBTK::PoissonSolver.

Member Data Documentation

◆ d_reinitializing_solver

bool IBTK::PETScNewtonKrylovSolver::d_reinitializing_solver = false
private

◆ d_petsc_x

Vec IBTK::PETScNewtonKrylovSolver::d_petsc_x = nullptr
private

◆ d_petsc_b

Vec IBTK::PETScNewtonKrylovSolver::d_petsc_b = nullptr
private

◆ d_petsc_r

Vec IBTK::PETScNewtonKrylovSolver::d_petsc_r = nullptr
private

◆ d_options_prefix

std::string IBTK::PETScNewtonKrylovSolver::d_options_prefix
private

◆ d_petsc_comm

MPI_Comm IBTK::PETScNewtonKrylovSolver::d_petsc_comm = PETSC_COMM_WORLD
private

◆ d_petsc_snes

SNES IBTK::PETScNewtonKrylovSolver::d_petsc_snes = nullptr
private

◆ d_petsc_jac

Mat IBTK::PETScNewtonKrylovSolver::d_petsc_jac = nullptr
private

◆ d_managing_petsc_snes

bool IBTK::PETScNewtonKrylovSolver::d_managing_petsc_snes = true
private

◆ d_user_provided_function

bool IBTK::PETScNewtonKrylovSolver::d_user_provided_function = false
private

◆ d_user_provided_jacobian

bool IBTK::PETScNewtonKrylovSolver::d_user_provided_jacobian = false
private

◆ d_F

SAMRAI::tbox::Pointer<GeneralOperator> IBTK::NewtonKrylovSolver::d_F
protectedinherited

◆ d_J

SAMRAI::tbox::Pointer<JacobianOperator> IBTK::NewtonKrylovSolver::d_J
protectedinherited

◆ d_krylov_solver

SAMRAI::tbox::Pointer<KrylovLinearSolver> IBTK::NewtonKrylovSolver::d_krylov_solver
protectedinherited

◆ d_x

SAMRAI::tbox::Pointer<SAMRAI::solv::SAMRAIVectorReal<NDIM, double> > IBTK::NewtonKrylovSolver::d_x
protectedinherited

◆ d_b

SAMRAI::tbox::Pointer<SAMRAI::solv::SAMRAIVectorReal<NDIM, double> > IBTK::NewtonKrylovSolver::d_b
protectedinherited

◆ d_r

SAMRAI::tbox::Pointer<SAMRAI::solv::SAMRAIVectorReal<NDIM, double> > IBTK::NewtonKrylovSolver::d_r
protectedinherited

◆ d_max_evaluations

int IBTK::NewtonKrylovSolver::d_max_evaluations = 10000
protectedinherited

◆ d_solution_tol

double IBTK::NewtonKrylovSolver::d_solution_tol = 1.0e-8
protectedinherited

◆ d_current_linear_iterations

int IBTK::NewtonKrylovSolver::d_current_linear_iterations = 0
protectedinherited

◆ d_object_name

std::string IBTK::GeneralSolver::d_object_name = "unitialized"
protectedinherited

◆ d_is_initialized

bool IBTK::GeneralSolver::d_is_initialized = false
protectedinherited

◆ d_homogeneous_bc

bool IBTK::GeneralSolver::d_homogeneous_bc = false
protectedinherited

◆ d_solution_time

double IBTK::GeneralSolver::d_solution_time = std::numeric_limits<double>::quiet_NaN()
protectedinherited

◆ d_current_time

double IBTK::GeneralSolver::d_current_time = std::numeric_limits<double>::quiet_NaN()
protectedinherited

◆ d_new_time

double IBTK::GeneralSolver::d_new_time = std::numeric_limits<double>::quiet_NaN()
protectedinherited

◆ d_rel_residual_tol

double IBTK::GeneralSolver::d_rel_residual_tol = 0.0
protectedinherited

◆ d_abs_residual_tol

double IBTK::GeneralSolver::d_abs_residual_tol = 0.0
protectedinherited

◆ d_max_iterations

int IBTK::GeneralSolver::d_max_iterations = 100
protectedinherited

◆ d_current_iterations

int IBTK::GeneralSolver::d_current_iterations = 0
protectedinherited

◆ d_current_residual_norm

double IBTK::GeneralSolver::d_current_residual_norm = std::numeric_limits<double>::quiet_NaN()
protectedinherited

◆ d_hier_math_ops

SAMRAI::tbox::Pointer<HierarchyMathOps> IBTK::GeneralSolver::d_hier_math_ops
protectedinherited

◆ d_hier_math_ops_external

bool IBTK::GeneralSolver::d_hier_math_ops_external = false
protectedinherited

◆ d_enable_logging

bool IBTK::GeneralSolver::d_enable_logging = false
protectedinherited

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