#include <source/solvers/packages/petsc/SNES_SAMRAIContext.h>
Inheritance diagram for SAMRAI::solv::SNES_SAMRAIContext< DIM >:
This class declares five private static member functions to link user-defined routines for nonlinear residual calculation, Jacobian evaluation, preconditioner setup and solve, and Jacobian-vector product operations. The implementation of these functions is defined by the user in a subclass of the abstract base class SNESAbstractFunctions. The vector objects used within the solver are provided by the PETSc_SAMRAIVectorReal<DIM> wrapper class.
If no parameters are read from input, PETSc defaults are used. See the PETSc documentation (http://www-unix.mcs.anl.gov/petsc/). for more information on default parameters and SNES functionality.
Optional input keys and types that can be read by this class are:
modifiedgramschmidt gmres_cgs_refine_ifneeded gmres_cgs_refine_always
See the PETSc documentation for more information.
Note that all input values may override values read in from restart. If no new input value is given, the restart value is used.
A sample input file entry might look like:
absolute_tolerance = 10.e-10
relative_tolerance = 10.e-6
step_tolerance = 10.e-8
maximum_nonlinear_iterations = 200
forcing_term_strategy = "EWCHOICE1"
Note that input values can also be set using accessor functions. Values that are set via this mechanism will be cached both in the solver context as well as in the corresponding PETSc object. Thus values changed on-the-fly will be written to restart. These input values can also be changed by directly accessing the corresponding PETSc object and using native PETSc function calls; however such settings/changes will NOT be cached in the solver context, and so will not be written to restart.
SAMRAI::solv::SNES_SAMRAIContext< DIM >::SNES_SAMRAIContext | ( | const std::string & | object_name, | |
tbox::Pointer< tbox::Database > | input_db, | |||
SNESAbstractFunctions * | my_functions | |||
) |
Constructor for SNES_SAMRAIContext<DIM> allocates the SNES object and initializes rudimentary state associated with user-supplied solver components. Then, it reads solver parameter from input and restart which may override default values.
When assertion checking is active, an unrecoverable exception will result if the name string is empty or the pointer to the user-defined SNES functions object is null.
SAMRAI::solv::SNES_SAMRAIContext< DIM >::~SNES_SAMRAIContext | ( | ) |
Destructor for solve_SNES_SAMRAIContext destroys the SNES and the PETSc solution vector wrapper.
SNES SAMRAI::solv::SNES_SAMRAIContext< DIM >::getSNESSolver | ( | ) | const |
Return the PETSc nonlinear solver object.
SNESAbstractFunctions * SAMRAI::solv::SNES_SAMRAIContext< DIM >::getSNESFunctions | ( | ) | const |
Return pointer to object providing user-defined functions for SNES.
KSP SAMRAI::solv::SNES_SAMRAIContext< DIM >::getKrylovSolver | ( | ) | const |
Return the PETSc Krylov solver object.
Mat SAMRAI::solv::SNES_SAMRAIContext< DIM >::getJacobianMatrix | ( | ) | const |
Return the PETSc Mat object for the Jacobian.
double SAMRAI::solv::SNES_SAMRAIContext< DIM >::getAbsoluteTolerance | ( | ) | const |
Get absolute tolerance for nonlinear solver.
void SAMRAI::solv::SNES_SAMRAIContext< DIM >::setAbsoluteTolerance | ( | double | abs_tol | ) |
Set absolute tolerance for nonlinear solver.
double SAMRAI::solv::SNES_SAMRAIContext< DIM >::getRelativeTolerance | ( | ) | const |
Get relative tolerance for nonlinear solver.
void SAMRAI::solv::SNES_SAMRAIContext< DIM >::setRelativeTolerance | ( | double | rel_tol | ) |
Set relative tolerance for nonlinear solver.
double SAMRAI::solv::SNES_SAMRAIContext< DIM >::getStepTolerance | ( | ) | const |
Get step tolerance for nonlinear solver.
void SAMRAI::solv::SNES_SAMRAIContext< DIM >::setStepTolerance | ( | double | step_tol | ) |
Set step tolerance for nonlinear solver.
int SAMRAI::solv::SNES_SAMRAIContext< DIM >::getMaxNonlinearIterations | ( | ) | const |
Get maximum iterations for nonlinear solver.
void SAMRAI::solv::SNES_SAMRAIContext< DIM >::setMaxNonlinearIterations | ( | int | max_nli | ) |
Set maximum iterations for nonlinear solver.
int SAMRAI::solv::SNES_SAMRAIContext< DIM >::getMaxFunctionEvaluations | ( | ) | const |
Get maximum function evaluations by nonlinear solver.
void SAMRAI::solv::SNES_SAMRAIContext< DIM >::setMaxFunctionEvaluations | ( | int | max_feval | ) |
Set maximum function evaluations in nonlinear solver.
std::string SAMRAI::solv::SNES_SAMRAIContext< DIM >::getForcingTermStrategy | ( | ) | const |
Get strategy for forcing term.
void SAMRAI::solv::SNES_SAMRAIContext< DIM >::setForcingTermStrategy | ( | std::string & | strategy | ) |
Set strategy for forcing term.
double SAMRAI::solv::SNES_SAMRAIContext< DIM >::getConstantForcingTerm | ( | ) | const |
Get value of constant forcing term.
void SAMRAI::solv::SNES_SAMRAIContext< DIM >::setConstantForcingTerm | ( | double | fixed_eta | ) |
Set value of constant forcing term.
double SAMRAI::solv::SNES_SAMRAIContext< DIM >::getInitialForcingTerm | ( | ) | const |
Get value of initial forcing term.
void SAMRAI::solv::SNES_SAMRAIContext< DIM >::setInitialForcingTerm | ( | double | initial_eta | ) |
Set value of initial forcing term.
double SAMRAI::solv::SNES_SAMRAIContext< DIM >::getMaximumForcingTerm | ( | ) | const |
Get value of maximum forcing term.
void SAMRAI::solv::SNES_SAMRAIContext< DIM >::setMaximumForcingTerm | ( | double | max_eta | ) |
Set value of maximum forcing term.
double SAMRAI::solv::SNES_SAMRAIContext< DIM >::getEWChoice2Exponent | ( | ) | const |
Get value of exponent in Eisenstat-Walker Choice 2 forcing term.
void SAMRAI::solv::SNES_SAMRAIContext< DIM >::setEWChoice2Exponent | ( | double | alpha | ) |
Set value of exponent in Eisenstat-Walker Choice 2 forcing term.
double SAMRAI::solv::SNES_SAMRAIContext< DIM >::getEWChoice2SafeguardExponent | ( | ) | const |
Get value of exponent in Eisenstat-Walker Choice 2 safeguard.
void SAMRAI::solv::SNES_SAMRAIContext< DIM >::setEWChoice2SafeguardExponent | ( | double | beta | ) |
Set value of exponent in Eisenstat-Walker Choice 2 safeguard.
double SAMRAI::solv::SNES_SAMRAIContext< DIM >::getEWChoice2ScaleFactor | ( | ) | const |
Get value of factor used to scale Eisenstat-Walker Choice 2 forcing term.
void SAMRAI::solv::SNES_SAMRAIContext< DIM >::setEWChoice2ScaleFactor | ( | double | gamma | ) |
Set value of factor used to scale Eisenstat-Walker Choice 2 forcing term.
double SAMRAI::solv::SNES_SAMRAIContext< DIM >::getEWSafeguardThreshold | ( | ) | const |
Get value of threshold to disable safeguard in Eisenstat-Walker forcing terms.
void SAMRAI::solv::SNES_SAMRAIContext< DIM >::setEWSafeguardThreshold | ( | double | threshold | ) |
Set value of threshold to disable safeguard in Eisenstat-Walker forcing terms.
std::string SAMRAI::solv::SNES_SAMRAIContext< DIM >::getLinearSolverType | ( | ) | const |
Get type of linear solver.
void SAMRAI::solv::SNES_SAMRAIContext< DIM >::setLinearSolverType | ( | std::string & | type | ) |
Set type of linear solver.
bool SAMRAI::solv::SNES_SAMRAIContext< DIM >::getUsesPreconditioner | ( | ) | const |
Get whether a preconditioner is used.
void SAMRAI::solv::SNES_SAMRAIContext< DIM >::setUsesPreconditioner | ( | bool | uses_preconditioner | ) |
Set whether to use a preconditioner.
double SAMRAI::solv::SNES_SAMRAIContext< DIM >::getLinearSolverAbsoluteTolerance | ( | ) | const |
Get absolute tolerance for linear solver.
void SAMRAI::solv::SNES_SAMRAIContext< DIM >::setLinearSolverAbsoluteTolerance | ( | double | abs_tol | ) |
Set absolute tolerance for linear solver.
double SAMRAI::solv::SNES_SAMRAIContext< DIM >::getLinearSolverDivergenceTolerance | ( | ) | const |
Get divergence tolerance for linear solver.
void SAMRAI::solv::SNES_SAMRAIContext< DIM >::setLinearSolverDivergenceTolerance | ( | double | div_tol | ) |
Set divergence tolerance for linear solver.
int SAMRAI::solv::SNES_SAMRAIContext< DIM >::getMaximumLinearIterations | ( | ) | const |
Get maximum linear iterations for linear solver.
void SAMRAI::solv::SNES_SAMRAIContext< DIM >::setMaximumLinearIterations | ( | int | max_li | ) |
Set maximum linear iterations for linear solver.
int SAMRAI::solv::SNES_SAMRAIContext< DIM >::getMaximumGMRESKrylovDimension | ( | ) | const |
Get maximum Krylov dimension in GMRES linear solver.
void SAMRAI::solv::SNES_SAMRAIContext< DIM >::setMaximumGMRESKrylovDimension | ( | int | d | ) |
Set maximum Krylov dimension in GMRES linear solver.
std::string SAMRAI::solv::SNES_SAMRAIContext< DIM >::getGMRESOrthogonalizationMethod | ( | ) | const |
Get orthogonalization method used in GMRES linear solver.
void SAMRAI::solv::SNES_SAMRAIContext< DIM >::setGMRESOrthogonalizationMethod | ( | std::string & | method | ) |
Set orthogonalization method used in GMRES linear solver.
bool SAMRAI::solv::SNES_SAMRAIContext< DIM >::getUsesExplicitJacobian | ( | ) | const |
Get whether a method for explicit Jacobian-vector products is provided.
void SAMRAI::solv::SNES_SAMRAIContext< DIM >::setUsesExplicitJacobian | ( | bool | use_jac | ) |
Set whether a method for explicit Jacobian-vector products is provided.
std::string SAMRAI::solv::SNES_SAMRAIContext< DIM >::getDifferencingParameterMethod | ( | ) | const |
Get method for computing differencing parameter.
void SAMRAI::solv::SNES_SAMRAIContext< DIM >::setDifferencingParameterMethod | ( | std::string & | method | ) |
Set method for computing differencing parameter.
double SAMRAI::solv::SNES_SAMRAIContext< DIM >::getFunctionEvaluationError | ( | ) | const |
Get estimate of error in function evaluation.
void SAMRAI::solv::SNES_SAMRAIContext< DIM >::setFunctionEvaluationError | ( | double | evaluation_error | ) |
Set estimate of error in function evaluation.
void SAMRAI::solv::SNES_SAMRAIContext< DIM >::initialize | ( | tbox::Pointer< SAMRAIVectorReal< DIM, double > > | solution | ) | [virtual] |
Initialize the state of the SNES solver based on vector argument representing the solution of the nonlinear system. In general, this routine must be called before the solve() routine is invoked.
Implements SAMRAI::solv::NonlinearSolverStrategy< DIM >.
void SAMRAI::solv::SNES_SAMRAIContext< DIM >::resetSolver | ( | const int | coarsest_level, | |
const int | finest_level | |||
) |
Reset the state of the nonlinear solver after regridding.
int SAMRAI::solv::SNES_SAMRAIContext< DIM >::solve | ( | ) | [virtual] |
Solve the nonlinear problem. In general, the initialize() routine must be called before this solve function to set up the solver. Returns 1 if successful, 0 otherwise.
Implements SAMRAI::solv::NonlinearSolverStrategy< DIM >.
int SAMRAI::solv::SNES_SAMRAIContext< DIM >::getNonlinearIterationCount | ( | ) | const |
Obtain number of nonlinear iterations.
int SAMRAI::solv::SNES_SAMRAIContext< DIM >::getTotalLinearIterationCount | ( | ) | const |
Obtain total number of linear iterations accumulated over all nonlinear iterations.
void SAMRAI::solv::SNES_SAMRAIContext< DIM >::reportCompletionCode | ( | std::ostream & | os = tbox::plog |
) | const |
Report reason for termination.
void SAMRAI::solv::SNES_SAMRAIContext< DIM >::putToDatabase | ( | tbox::Pointer< tbox::Database > | db | ) | [virtual] |
Write solver parameters to restart database matching object name.
When assertion checking is active, an unrecoverable exception will result if database pointer is null.
Implements SAMRAI::tbox::Serializable.
void SAMRAI::solv::SNES_SAMRAIContext< DIM >::printClassData | ( | std::ostream & | os | ) | const [virtual] |
Print out all members of integrator instance to given output stream.