|
IBAMR
IBAMR version 0.19.
|
Class PETScKrylovLinearSolver provides a KrylovLinearSolver interface for a PETSc Krylov subspace iterative linear solver (KSP). More...
#include <ibtk/PETScKrylovLinearSolver.h>

Public Member Functions | |
| PETScKrylovLinearSolver (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 KrylovLinearSolver that employs the PETSc KSP solver framework. More... | |
| PETScKrylovLinearSolver (std::string object_name, const KSP &petsc_ksp) | |
| Constructor for a concrete KrylovLinearSolver that acts as a "wrapper" for a supplied PETSc KSP object. More... | |
| ~PETScKrylovLinearSolver () | |
| Destructor. More... | |
| void | setKSPType (const std::string &ksp_type) |
| Set the KSP type. 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... | |
Static Public Member Functions | |
| static SAMRAI::tbox::Pointer< KrylovLinearSolver > | 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 PETScKrylovLinearSolver. More... | |
Friends | |
| class | PETScNewtonKrylovSolver |
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... | |
| 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) |
Krylov solver functionality. | |
| virtual SAMRAI::tbox::Pointer< LinearOperator > | getOperator () const |
| Retrieve the linear operator used when solving \(Ax=b\). More... | |
| virtual SAMRAI::tbox::Pointer< LinearSolver > | getPreconditioner () const |
| Retrieve the preconditioner used by the Krylov subspace method when solving \(Ax=b\). More... | |
| SAMRAI::tbox::Pointer< LinearOperator > | d_A |
| SAMRAI::tbox::Pointer< LinearSolver > | d_pc_solver |
| SAMRAI::tbox::Pointer< SAMRAI::solv::SAMRAIVectorReal< NDIM, double > > | d_x |
| SAMRAI::tbox::Pointer< SAMRAI::solv::SAMRAIVectorReal< NDIM, double > > | d_b |
Functions to access solver parameters. | |
| virtual void | setInitialGuessNonzero (bool initial_guess_nonzero=true) |
| Set whether the initial guess is non-zero. More... | |
| virtual bool | getInitialGuessNonzero () const |
| Get whether the initial guess is non-zero. More... | |
Logging functions. | |
| virtual void | printClassData (std::ostream &stream) override |
| Print class data to stream. More... | |
| bool | d_initial_guess_nonzero = true |
| bool | d_nullspace_contains_constant_vec = false |
| std::vector< SAMRAI::tbox::Pointer< SAMRAI::solv::SAMRAIVectorReal< NDIM, double > > > | d_nullspace_basis_vecs |
Krylov solver functionality. | |
| void | setOperator (SAMRAI::tbox::Pointer< LinearOperator > A) override |
| Set the linear operator used when solving \(Ax=b\). More... | |
| void | setPreconditioner (SAMRAI::tbox::Pointer< LinearSolver > pc_solver=nullptr) override |
| Set the preconditioner used by the Krylov subspace method when solving \(Ax=b\). More... | |
| void | setNullSpace (bool contains_constant_vec, const std::vector< SAMRAI::tbox::Pointer< SAMRAI::solv::SAMRAIVectorReal< NDIM, double > > > &nullspace_basis_vecs=std::vector< SAMRAI::tbox::Pointer< SAMRAI::solv::SAMRAIVectorReal< NDIM, double > > >()) override |
| Set the nullspace of the linear system. More... | |
| bool | solveSystem (SAMRAI::solv::SAMRAIVectorReal< NDIM, double > &x, SAMRAI::solv::SAMRAIVectorReal< NDIM, double > &b) override |
| Solve the linear system of equations \(Ax=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 \(Ax=b\). More... | |
| void | deallocateSolverState () override |
| Remove all hierarchy dependent data allocated by initializeSolverState(). More... | |
| PETScKrylovLinearSolver ()=delete | |
| Default constructor. More... | |
| PETScKrylovLinearSolver (const PETScKrylovLinearSolver &from)=delete | |
| Copy constructor. More... | |
| PETScKrylovLinearSolver & | operator= (const PETScKrylovLinearSolver &that)=delete |
| Assignment operator. More... | |
| void | common_ctor () |
| Common routine used by all class constructors. More... | |
| void | resetWrappedKSP (KSP &petsc_ksp) |
| Reset the KSP wrapped by this solver class. More... | |
| void | resetKSPOptions () |
| Reset the values of the convergence tolerances for the PETSc KSP object. More... | |
| void | resetKSPOperators () |
| Reset the KSP operators to correspond to the supplied LinearOperator. More... | |
| void | resetKSPPC () |
| Reset the KSP PC to correspond to the supplied preconditioner. More... | |
| void | resetMatNullSpace () |
| Reset the Mat nullspace object to correspond to the supplied nullspace basis vectors. More... | |
| void | deallocateNullSpaceData () |
| Destroy data allocated to describe nullspace. More... | |
Static functions for use by PETSc KSP and MatShell objects. | |
| std::string | d_ksp_type |
| bool | d_reinitializing_solver = false |
| Vec | d_petsc_x = nullptr |
| Vec | d_petsc_b = nullptr |
| std::string | d_options_prefix |
| MPI_Comm | d_petsc_comm |
| KSP | d_petsc_ksp = nullptr |
| Mat | d_petsc_mat = nullptr |
| MatNullSpace | d_petsc_nullsp = nullptr |
| bool | d_managing_petsc_ksp = true |
| bool | d_user_provided_mat = false |
| bool | d_user_provided_pc = false |
| SAMRAI::tbox::Pointer< SAMRAI::solv::SAMRAIVectorReal< NDIM, double > > | d_nullspace_constant_vec |
| Vec | d_petsc_nullspace_constant_vec = nullptr |
| std::vector< Vec > | d_petsc_nullspace_basis_vecs |
| bool | d_solver_has_attached_nullspace = false |
| static PetscErrorCode | MatVecMult_SAMRAI (Mat A, Vec x, Vec y) |
| Compute the matrix vector product \(y=Ax\). More... | |
| static PetscErrorCode | PCApply_SAMRAI (PC pc, Vec x, Vec y) |
| Apply the preconditioner to x and store the result in y. More... | |
This solver class provides access to a large number of Krylov subspace solvers for linear problems of the form \(Ax=b\) using the PETSc KSP linear solver interface. See http://www.mcs.anl.gov/petsc/documentation/linearsolvertable.html for a complete list of the Krylov solvers provided by this class. 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() ksp_type = "gmres" // see setKSPType() initial_guess_nonzero = TRUE // see setInitialGuessNonzero() rel_residual_tol = 1.0e-5 // see setRelativeTolerance() abs_residual_tol = 1.0e-50 // see setAbsoluteTolerance() max_iterations = 10000 // 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::PETScKrylovLinearSolver::PETScKrylovLinearSolver | ( | 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::PETScKrylovLinearSolver::PETScKrylovLinearSolver | ( | std::string | object_name, |
| const KSP & | petsc_ksp | ||
| ) |
| object_name | Name of the solver |
| petsc_ksp | PETSc KSP object |
| IBTK::PETScKrylovLinearSolver::~PETScKrylovLinearSolver | ( | ) |
|
privatedelete |
|
privatedelete |
| from | The value to copy to this object. |
|
inlinestatic |
| void IBTK::PETScKrylovLinearSolver::setKSPType | ( | const std::string & | ksp_type | ) |
| void IBTK::PETScKrylovLinearSolver::setOptionsPrefix | ( | const std::string & | options_prefix | ) |
| const KSP& IBTK::PETScKrylovLinearSolver::getPETScKSP | ( | ) | const |
|
overridevirtual |
Reimplemented from IBTK::KrylovLinearSolver.
|
overridevirtual |
Reimplemented from IBTK::KrylovLinearSolver.
|
overridevirtual |
Basis vectors must be orthogonal but are not required to be orthonormal. Basis vectors will be normalized automatically.
Reimplemented from IBTK::LinearSolver.
|
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 linear operator or preconditioner objects have been registered with this class via setOperator() and setPreconditioner(), 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 linear operator or preconditioner objects have been registered with this class via setOperator() and setPreconditioner(), 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 |
|
private |
|
private |
|
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 |
Reimplemented in IBTK::PETScPCLSWrapper.
|
overridevirtualinherited |
Reimplemented from IBTK::GeneralSolver.
|
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 |
|
protectedinherited |
|
protectedvirtualinherited |
Reimplemented in IBTK::PoissonSolver.
|
friend |
|
private |
|
private |
|
private |
|
private |
|
private |
|
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 |
1.8.17