Reference documentation for deal.II version 9.5.0
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
Public Types | Public Member Functions | Public Attributes | List of all members
SUNDIALS::KINSOL< VectorType >::AdditionalData Class Reference

#include <deal.II/sundials/kinsol.h>

Public Types

enum  SolutionStrategy { newton = KIN_NONE , linesearch = KIN_LINESEARCH , fixed_point = KIN_FP , picard = KIN_PICARD }
 

Public Member Functions

 AdditionalData (const SolutionStrategy &strategy=linesearch, const unsigned int maximum_non_linear_iterations=200, const double function_tolerance=0.0, const double step_tolerance=0.0, const bool no_init_setup=false, const unsigned int maximum_setup_calls=0, const double maximum_newton_step=0.0, const double dq_relative_error=0.0, const unsigned int maximum_beta_failures=0, const unsigned int anderson_subspace_size=0)
 
void add_parameters (ParameterHandler &prm)
 

Public Attributes

SolutionStrategy strategy
 
unsigned int maximum_non_linear_iterations
 
double function_tolerance
 
double step_tolerance
 
bool no_init_setup
 
unsigned int maximum_setup_calls
 
double maximum_newton_step
 
double dq_relative_error
 
unsigned int maximum_beta_failures
 
unsigned int anderson_subspace_size
 

Detailed Description

template<typename VectorType = Vector<double>>
class SUNDIALS::KINSOL< VectorType >::AdditionalData

Additional parameters that can be passed to the KINSOL class.

Definition at line 190 of file kinsol.h.

Member Enumeration Documentation

◆ SolutionStrategy

template<typename VectorType = Vector<double>>
enum SUNDIALS::KINSOL::AdditionalData::SolutionStrategy

KINSOL solution strategy. KINSOL includes a Newton-Krylov solver (both local and global) as well as Picard and fixed point solvers.

Enumerator
newton 

Standard Newton iteration.

linesearch 

Newton iteration with linesearch.

fixed_point 

Fixed point iteration.

picard 

Picard iteration.

Definition at line 197 of file kinsol.h.

Constructor & Destructor Documentation

◆ AdditionalData()

template<typename VectorType >
SUNDIALS::KINSOL< VectorType >::AdditionalData::AdditionalData ( const SolutionStrategy strategy = linesearch,
const unsigned int  maximum_non_linear_iterations = 200,
const double  function_tolerance = 0.0,
const double  step_tolerance = 0.0,
const bool  no_init_setup = false,
const unsigned int  maximum_setup_calls = 0,
const double  maximum_newton_step = 0.0,
const double  dq_relative_error = 0.0,
const unsigned int  maximum_beta_failures = 0,
const unsigned int  anderson_subspace_size = 0 
)

Initialization parameters for KINSOL.

Global parameters:

Parameters
strategySolution strategy
maximum_non_linear_iterationsMaximum number of nonlinear iterations
function_toleranceFunction norm stopping tolerance
step_toleranceScaled step stopping tolerance

Newton parameters:

Parameters
no_init_setupNo initial matrix setup
maximum_setup_callsMaximum iterations without matrix setup
maximum_newton_stepMaximum allowable scaled length of the Newton step
dq_relative_errorRelative error for different quotient computation

Line search parameters:

Parameters
maximum_beta_failuresMaximum number of beta-condition failures

Fixed point and Picard parameters:

Parameters
anderson_subspace_sizeAnderson acceleration subspace size

Definition at line 57 of file kinsol.cc.

Member Function Documentation

◆ add_parameters()

template<typename VectorType >
void SUNDIALS::KINSOL< VectorType >::AdditionalData::add_parameters ( ParameterHandler prm)

Add all AdditionalData() parameters to the given ParameterHandler object. When the parameters are parsed from a file, the internal parameters are automatically updated.

The following parameters are declared:

set Function norm stopping tolerance = 0
set Maximum number of nonlinear iterations = 200
set Scaled step stopping tolerance = 0
set Solution strategy = linesearch
subsection Fixed point and Picard parameters
set Anderson acceleration subspace size = 5
end
subsection Linesearch parameters
set Maximum number of beta-condition failures = 0
end
subsection Newton parameters
set Maximum allowable scaled length of the Newton step = 0
set Maximum iterations without matrix setup = 0
set No initial matrix setup = false
set Relative error for different quotient computation = 0
end
SolutionStrategy strategy
Definition kinsol.h:304

These are one-to-one with the options you can pass at construction time.

The options you pass at construction time are set as default values in the ParameterHandler object prm. You can later modify them by parsing a parameter file using prm. The values of the parameter will be updated whenever the content of prm is updated.

Make sure that this class lives longer than prm. Undefined behavior will occur if you destroy this class, and then parse a parameter file using prm.

Definition at line 84 of file kinsol.cc.

Member Data Documentation

◆ strategy

template<typename VectorType = Vector<double>>
SolutionStrategy SUNDIALS::KINSOL< VectorType >::AdditionalData::strategy

The solution strategy to use. If you choose SolutionStrategy::newton, SolutionStrategy::linesearch, or SolutionStrategy::picard you have to provide the function residual(). If you choose SolutionStrategy::fixed_point, you have to provide the function iteration_function().

Definition at line 304 of file kinsol.h.

◆ maximum_non_linear_iterations

template<typename VectorType = Vector<double>>
unsigned int SUNDIALS::KINSOL< VectorType >::AdditionalData::maximum_non_linear_iterations

Maximum number of nonlinear iterations allowed.

Definition at line 309 of file kinsol.h.

◆ function_tolerance

template<typename VectorType = Vector<double>>
double SUNDIALS::KINSOL< VectorType >::AdditionalData::function_tolerance

A scalar used as a stopping tolerance on the scaled maximum norm of the system function \(F(u)\) or \(G(u)\).

If set to zero, default values provided by KINSOL will be used.

Definition at line 317 of file kinsol.h.

◆ step_tolerance

template<typename VectorType = Vector<double>>
double SUNDIALS::KINSOL< VectorType >::AdditionalData::step_tolerance

A scalar used as a stopping tolerance on the minimum scaled step length.

If set to zero, default values provided by KINSOL will be used.

Definition at line 325 of file kinsol.h.

◆ no_init_setup

template<typename VectorType = Vector<double>>
bool SUNDIALS::KINSOL< VectorType >::AdditionalData::no_init_setup

Whether an initial call to the preconditioner or Jacobian setup function should be made or not.

A call to this function is useful when solving a sequence of problems, in which the final preconditioner or Jacobian value from one problem is to be used initially for the next problem.

Definition at line 335 of file kinsol.h.

◆ maximum_setup_calls

template<typename VectorType = Vector<double>>
unsigned int SUNDIALS::KINSOL< VectorType >::AdditionalData::maximum_setup_calls

The maximum number of nonlinear iterations that can be performed between calls to the setup_jacobian() function.

If set to zero, default values provided by KINSOL will be used, and in practice this often means that KINSOL will re-use a Jacobian matrix computed in one iteration for later iterations.

Definition at line 345 of file kinsol.h.

◆ maximum_newton_step

template<typename VectorType = Vector<double>>
double SUNDIALS::KINSOL< VectorType >::AdditionalData::maximum_newton_step

The maximum allowable scaled length of the Newton step.

If set to zero, default values provided by KINSOL will be used.

Definition at line 352 of file kinsol.h.

◆ dq_relative_error

template<typename VectorType = Vector<double>>
double SUNDIALS::KINSOL< VectorType >::AdditionalData::dq_relative_error

The relative error in computing \(F(u)\), which is used in the difference quotient approximation to the Jacobian matrix when the user does not supply a solve_with_jacobian() function.

If set to zero, default values provided by KINSOL will be used.

Definition at line 361 of file kinsol.h.

◆ maximum_beta_failures

template<typename VectorType = Vector<double>>
unsigned int SUNDIALS::KINSOL< VectorType >::AdditionalData::maximum_beta_failures

The maximum number of beta-condition failures in the linesearch algorithm. Only used if strategy==SolutionStrategy::linesearch.

Definition at line 368 of file kinsol.h.

◆ anderson_subspace_size

template<typename VectorType = Vector<double>>
unsigned int SUNDIALS::KINSOL< VectorType >::AdditionalData::anderson_subspace_size

The size of the subspace used with Anderson acceleration in conjunction with Picard or fixed-point iteration.

If you set this to 0, no acceleration is used.

Definition at line 376 of file kinsol.h.


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