#include <source/solvers/packages/sundials/kinsol/KINSOLSolver.h>
Public Member Functions | |
KINSOLSolver (const std::string &object_name, KINSOLAbstractFunctions *my_functions, const int uses_preconditioner, const int uses_jac_times_vector) | |
virtual | ~KINSOLSolver () |
void | initialize (SundialsAbstractVector *solution, SundialsAbstractVector *uscale=NULL, SundialsAbstractVector *fscale=NULL) |
int | solve () |
void | setLogFileData (const std::string &log_fname, const int flag) |
void | setKINSOLFunctions (KINSOLAbstractFunctions *my_functions, const int uses_preconditioner, const int uses_jac_times_vector) |
void | setPreconditioner (const int uses_preconditioner) |
void | setJacobianTimesVector (const int uses_jac_times_vector) |
KINSOLAbstractFunctions * | getKINSOLFunctions () const |
void | setConstraintVector (SundialsAbstractVector *constraints) |
void | setResidualStoppingTolerance (const double tol) |
void | setMaxIterations (const int maxits) |
void | setMaxKrylovDimension (const int kdim) |
void | setGlobalStrategy (const int global) |
void | setMaxNewtonStep (const double maxstep) |
void | setNonlinearStepTolerance (const double tol) |
void | setRelativeFunctionError (const double reserr) |
void | setLinearSolverConvergenceTest (const int conv) |
void | setLinearSolverConstantTolerance (const double tol) |
void | setEisenstatWalkerParameters (const double alpha, const double gamma) |
void | setMaxStepsWithNoPrecondSetup (const int maxsolv) |
void | setMaxLinearSolveRestarts (const int restarts) |
void | setMaxSubSetupCalls (const int maxsub) |
void | setResidualMonitoringParams (const double omega_min, const double omega_max) |
void | setResidualMonitoringConstant (const double omega) |
void | setNoMinEps (const bool flag) |
void | setMaxBetaFails (const int max_beta_fails) |
void | setNoInitialSetup (const bool flag) |
void | setNoResidualMonitoring (const bool flag) |
int | getTotalNumberOfNonlinearIterations () const |
int | getTotalNumberOfFunctionCalls () const |
int | getTotalNumberOfBetaConditionFailures () const |
int | getTotalNumberOfBacktracks () const |
double | getScaledResidualNorm () const |
double | getNewtonStepLength () const |
virtual void | printClassData (std::ostream &os) const |
void | initializeKINSOL () |
Note that this class provides no input or restart capabilities and relies on KINSOL for output reporting. When using KINSOL in an application using SAMRAI, it is straightforward to include this functionality in the entity using this solver class.
KINSOL was developed in the Center for Applied Scientific Computing (CASC) at Lawrence Livermore National Laboratory (LLNL). For more information about KINSOL and a complete description of the operations and data structures used by this class, see A.G. Taylor and A.C. Hindmarsh, "User documentation for KINSOL, a nonlinear solver for sequential and parallel computers", UCRL-ID-131185, Lawrence Livermore National Laboratory, 1998.
SAMRAI::solv::KINSOLSolver::KINSOLSolver | ( | const std::string & | object_name, | |
KINSOLAbstractFunctions * | my_functions, | |||
const int | uses_preconditioner, | |||
const int | uses_jac_times_vector | |||
) |
Constructor for KINSOLSolver sets default KINSOL parameters and initializes the solver package with user-supplied functions. Solver parameters may be changed later using member functions described below. The integer flags indicate whether user-supplied preconditioner and Jacobian-vector product function should be used. Zero indicates no user function; otherwise, user function will be used by the nonlinear solver.
Important note: The solution vector is not passed into the constructor. Before the solver can be used, the initialize() function must be called.
When assertion checking is active, an unrecoverable assertion will result if pointer to functions is null or string is empty.
SAMRAI::solv::KINSOLSolver::~KINSOLSolver | ( | ) | [virtual] |
Virtual destructor for KINSOLSolver.
void SAMRAI::solv::KINSOLSolver::initialize | ( | SundialsAbstractVector * | solution, | |
SundialsAbstractVector * | uscale = NULL , |
|||
SundialsAbstractVector * | fscale = NULL | |||
) |
Initialize solver with solution vector. The solution vector is required to initialize the memory record used internally within KINSOL. This routine must be called before the solver can be used.
When assertion checking is active, an unrecoverable assertion will result if vector pointer is null.
Optionally set the scaling vectors used by KINSOL to scale either nonlinear solution vector or nonlinear residual vector. The elements of the scaling vectors must be positive. In either case, the scaling vector should be defined so that the vector formed by taking the element-wise product of the solution/residual vector and scaling vector has all elements roughly the same magnitude when the solution vector IS/IS NOT NEAR a root of the nonlinear function.
See KINSOL documentation for more information.
int SAMRAI::solv::KINSOLSolver::solve | ( | ) |
Solve nonlinear problem and return integer termination code defined by KINSOL. The default return value is KINSOL_SUCCESS (= 1) indicating success. Return values which indicate non-recoverable nonlinear solver behavior are KINSOL_NO_MEM (= -1), KINSOL_INPUT_ERROR (= -2), and KINSOL_LSOLV_NO_MEM (= -3). Return values PRECONDSET_FAILURE (= 9), and PRECONDSOLVE_FAILURE (= 10) generally indicate non-recoverable behavior in the preconditioner. See kinsol.h header file for more information about return values.
If KINSOL requires re-initialization, it is automatically done before the solve. This may be required if any of the KINSOL data parameters have changed since the last call to the solver.
void SAMRAI::solv::KINSOLSolver::setLogFileData | ( | const std::string & | log_fname, | |
const int | flag | |||
) |
Accessory function for setting KINSOL output log file name and output printing options. Output file name and options may be changed throughout run as desired.
KINSOL printing options are:
The default is no output (i.e., 0). If the file name string is empty the default file name "kinsol.log" is used.
See KINSOL documentation for more information.
void SAMRAI::solv::KINSOLSolver::setKINSOLFunctions | ( | KINSOLAbstractFunctions * | my_functions, | |
const int | uses_preconditioner, | |||
const int | uses_jac_times_vector | |||
) |
Accessory functions for passing user-defined function information to KINSOL.
my_functions is a pointer to the abstract function subclass object that defines the residual calculation and preconditioner functions.
uses_preconditioner turns user preconditioner on or off.
uses_jac_times_vector turns user Jacobian-vector product on or off.
Flags use "TRUE"/"FALSE" values defined in KINSOL. See KINSOL documentation for more information.
void SAMRAI::solv::KINSOLSolver::setPreconditioner | ( | const int | uses_preconditioner | ) |
void SAMRAI::solv::KINSOLSolver::setJacobianTimesVector | ( | const int | uses_jac_times_vector | ) |
KINSOLAbstractFunctions * SAMRAI::solv::KINSOLSolver::getKINSOLFunctions | ( | ) | const |
Return pointer to object that provides user-defined functions for KINSOL.
void SAMRAI::solv::KINSOLSolver::setConstraintVector | ( | SundialsAbstractVector * | constraints | ) |
Set constraints on nonlinear solution. By default the constraint vector is null.
The constraints are applied in KINSOL as follows:
See KINSOL documentation for more information.
void SAMRAI::solv::KINSOLSolver::setResidualStoppingTolerance | ( | const double | tol | ) |
Accessory functions for setting nonlinear solver parameters. Parameters and default values are:
Residual stopping tolerance is tolerarnce on max_norm(fscale * residual), where product of vectors is another vector each element of which is the product of the corresponding entries in the original vectors. The default is .
Default maximum nonlinear iterations is 200.
Default maximum Krylov dimension is 1.
Options for global Newton method are: INEXACT_NEWTON = 0, LINESEARCH = 1. The default is INEXACT_NEWTON.
Default maximum Newton step is 1000*max(norm(uscale*u_0), norm(uscale)), where u_0 is the initial guess at the solution.
Default scaled step tolerarnce between successive nonlinear iterates is .
Default relative error for nonlinear function is set to machine_epsilon.
Scalar update constraint value restricts update of solution to del(u)/u < constraint_value. Here, vector ratio is another vector each element of which is the ratio of the corresponding entries in the original vectors. The default is no constraint.
See KINSOL documentation for more information.
void SAMRAI::solv::KINSOLSolver::setMaxIterations | ( | const int | maxits | ) |
void SAMRAI::solv::KINSOLSolver::setMaxKrylovDimension | ( | const int | kdim | ) |
void SAMRAI::solv::KINSOLSolver::setGlobalStrategy | ( | const int | global | ) |
void SAMRAI::solv::KINSOLSolver::setMaxNewtonStep | ( | const double | maxstep | ) |
void SAMRAI::solv::KINSOLSolver::setNonlinearStepTolerance | ( | const double | tol | ) |
void SAMRAI::solv::KINSOLSolver::setRelativeFunctionError | ( | const double | reserr | ) |
void SAMRAI::solv::KINSOLSolver::setLinearSolverConvergenceTest | ( | const int | conv | ) |
Accessory functions for setting convergence tests for inner linear solvers within an inexact Newton method. In general, the linear solver attempts to produce a step p, satisfying: norm(F(u) + J(u)*p) <= (eta + u_round)*norm(F(u)), where the norm is a scaled L2-norm.
The convergence test indicates the value for eta; options are:
The default option is ETACONSTANT.
The default constant value for eta is 0.1.
For choice ETACHOICE2, alpha = 2.0 and gamma = 0.9 are defaults.
See KINSOL documentation for more information.
void SAMRAI::solv::KINSOLSolver::setLinearSolverConstantTolerance | ( | const double | tol | ) |
void SAMRAI::solv::KINSOLSolver::setEisenstatWalkerParameters | ( | const double | alpha, | |
const double | gamma | |||
) |
void SAMRAI::solv::KINSOLSolver::setMaxStepsWithNoPrecondSetup | ( | const int | maxsolv | ) |
void SAMRAI::solv::KINSOLSolver::setMaxLinearSolveRestarts | ( | const int | restarts | ) |
void SAMRAI::solv::KINSOLSolver::setMaxSubSetupCalls | ( | const int | maxsub | ) |
The number of nonlinear iterations between checks by the nonlinear residual monitoring algorithm (specifies lenght of subinterval) NOTE: should be a multiple of MaxStepsWithNoPrecondSetup
void SAMRAI::solv::KINSOLSolver::setResidualMonitoringParams | ( | const double | omega_min, | |
const double | omega_max | |||
) |
Set values of omega_min and omega_max scalars used by nonlinear residual monitoring algorithm.
Defaults is [0.00001 and 0.9]
void SAMRAI::solv::KINSOLSolver::setResidualMonitoringConstant | ( | const double | omega | ) |
Set constant value used by residual monitoring algorithm. If omega=0, then it is estimated using omega_min and omega_max.
Default is 0.0.
void SAMRAI::solv::KINSOLSolver::setNoMinEps | ( | const bool | flag | ) |
Set flag controlling whether or not the value * of eps is bounded below by 0.01*fnormtol.
FALSE constrain value of eps by setting to the following: eps = MAX{0.01*fnormtol, eps}
TRUE do notconstrain value of eps
Default is FALSE
void SAMRAI::solv::KINSOLSolver::setMaxBetaFails | ( | const int | max_beta_fails | ) |
Set maximum number of beta condition failures in the line search algorithm.
Default is [MXNBCF_DEFAULT] (defined in kinsol_impl.h)
void SAMRAI::solv::KINSOLSolver::setNoInitialSetup | ( | const bool | flag | ) |
Flag controlling whether or not the KINSol routine makes an initial call to the linearl solver setup routine. Default is false.
void SAMRAI::solv::KINSOLSolver::setNoResidualMonitoring | ( | const bool | flag | ) |
Flag controlling whether or not the nonlinear residual monitoring schemes is used to control Jacobian updating Default is FALSE.
int SAMRAI::solv::KINSOLSolver::getTotalNumberOfNonlinearIterations | ( | ) | const |
Accessory functions to retrieve information fom KINSOL.
See KINSOL documentation for more information.
int SAMRAI::solv::KINSOLSolver::getTotalNumberOfFunctionCalls | ( | ) | const |
int SAMRAI::solv::KINSOLSolver::getTotalNumberOfBetaConditionFailures | ( | ) | const |
int SAMRAI::solv::KINSOLSolver::getTotalNumberOfBacktracks | ( | ) | const |
double SAMRAI::solv::KINSOLSolver::getScaledResidualNorm | ( | ) | const |
double SAMRAI::solv::KINSOLSolver::getNewtonStepLength | ( | ) | const |
void SAMRAI::solv::KINSOLSolver::printClassData | ( | std::ostream & | os | ) | const [virtual] |
Print out all data members for this object.
void SAMRAI::solv::KINSOLSolver::initializeKINSOL | ( | ) |