Reference documentation for deal.II version GIT relicensing-214-g6e74dec06b 2024-03-27 18:10:01+00:00
\(\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::IDA< VectorType >::AdditionalData Class Reference

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

Public Types

enum  InitialConditionCorrection { none = 0 , use_y_diff = 1 , use_y_dot = 2 }
 

Public Member Functions

 AdditionalData (const double initial_time=0.0, const double final_time=1.0, const double initial_step_size=1e-2, const double output_period=1e-1, const double minimum_step_size=1e-6, const unsigned int maximum_order=5, const unsigned int maximum_non_linear_iterations=10, const double ls_norm_factor=0, const double absolute_tolerance=1e-6, const double relative_tolerance=1e-5, const bool ignore_algebraic_terms_for_errors=true, const InitialConditionCorrection &ic_type=use_y_diff, const InitialConditionCorrection &reset_type=use_y_diff, const unsigned int maximum_non_linear_iterations_ic=5)
 
void add_parameters (ParameterHandler &prm)
 

Public Attributes

double initial_time
 
double final_time
 
double initial_step_size
 
double minimum_step_size
 
double absolute_tolerance
 
double relative_tolerance
 
unsigned int maximum_order
 
double output_period
 
bool ignore_algebraic_terms_for_errors
 
InitialConditionCorrection ic_type
 
InitialConditionCorrection reset_type
 
unsigned maximum_non_linear_iterations_ic
 
unsigned int maximum_non_linear_iterations
 
double ls_norm_factor
 

Detailed Description

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

Additional parameters that can be passed to the IDA class.

Definition at line 234 of file ida.h.

Member Enumeration Documentation

◆ InitialConditionCorrection

template<typename VectorType = Vector<double>>
enum SUNDIALS::IDA::AdditionalData::InitialConditionCorrection

IDA is a Differential Algebraic solver. As such, it requires initial conditions also for the first order derivatives. If you do not provide consistent initial conditions, (i.e., conditions for which F(y_dot(0), y(0), 0) = 0), you can ask SUNDIALS to compute initial conditions for you by specifying InitialConditionCorrection for the initial conditions both at the initial_time (ic_type) and after a reset has occurred (reset_type).

Enumerator
none 

Do not try to make initial conditions consistent.

use_y_diff 

Compute the algebraic components of y and differential components of y_dot, given the differential components of y. This option requires that the user specifies differential and algebraic components in the function get_differential_components.

use_y_dot 

Compute all components of y, given y_dot.

Definition at line 246 of file ida.h.

Constructor & Destructor Documentation

◆ AdditionalData()

template<typename VectorType = Vector<double>>
SUNDIALS::IDA< VectorType >::AdditionalData::AdditionalData ( const double  initial_time = 0.0,
const double  final_time = 1.0,
const double  initial_step_size = 1e-2,
const double  output_period = 1e-1,
const double  minimum_step_size = 1e-6,
const unsigned int  maximum_order = 5,
const unsigned int  maximum_non_linear_iterations = 10,
const double  ls_norm_factor = 0,
const double  absolute_tolerance = 1e-6,
const double  relative_tolerance = 1e-5,
const bool  ignore_algebraic_terms_for_errors = true,
const InitialConditionCorrection ic_type = use_y_diff,
const InitialConditionCorrection reset_type = use_y_diff,
const unsigned int  maximum_non_linear_iterations_ic = 5 
)
inline

Initialization parameters for IDA.

Global parameters:

Parameters
initial_timeInitial time
final_timeFinal time
initial_step_sizeInitial step size
output_periodTime interval between each output

Running parameters:

Parameters
minimum_step_sizeMinimum step size
maximum_orderMaximum BDF order
maximum_non_linear_iterationsMaximum number of nonlinear iterations
ls_norm_factorConverting factor from the integrator tolerance to the linear solver tolerance iterations

Error parameters:

Parameters
absolute_toleranceAbsolute error tolerance
relative_toleranceRelative error tolerance
ignore_algebraic_terms_for_errorsIgnore algebraic terms for error computations

Initial condition correction parameters:

Parameters
ic_typeInitial condition correction type
reset_typeInitial condition correction type after restart
maximum_non_linear_iterations_icInitial condition Newton max iterations

Definition at line 301 of file ida.h.

Member Function Documentation

◆ add_parameters()

template<typename VectorType = Vector<double>>
void SUNDIALS::IDA< VectorType >::AdditionalData::add_parameters ( ParameterHandler prm)
inline

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 Final time = 1.000000
set Initial time = 0.000000
set Time interval between each output = 0.2
subsection Error control
set Absolute error tolerance = 0.000001
set Ignore algebraic terms for error computations = true
set Relative error tolerance = 0.00001
set Use local tolerances = false
end
subsection Initial condition correction parameters
set Correction type at initial time = none
set Correction type after restart = none
set Maximum number of nonlinear iterations = 5
end
subsection Running parameters
set Initial step size = 0.1
set Maximum number of nonlinear iterations = 10
set Maximum order of BDF = 5
set Minimum step size = 0.000001
end

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 378 of file ida.h.

Member Data Documentation

◆ initial_time

template<typename VectorType = Vector<double>>
double SUNDIALS::IDA< VectorType >::AdditionalData::initial_time

Initial time for the DAE.

Definition at line 463 of file ida.h.

◆ final_time

template<typename VectorType = Vector<double>>
double SUNDIALS::IDA< VectorType >::AdditionalData::final_time

Final time.

Definition at line 468 of file ida.h.

◆ initial_step_size

template<typename VectorType = Vector<double>>
double SUNDIALS::IDA< VectorType >::AdditionalData::initial_step_size

Initial step size.

Definition at line 473 of file ida.h.

◆ minimum_step_size

template<typename VectorType = Vector<double>>
double SUNDIALS::IDA< VectorType >::AdditionalData::minimum_step_size

Minimum step size.

Definition at line 478 of file ida.h.

◆ absolute_tolerance

template<typename VectorType = Vector<double>>
double SUNDIALS::IDA< VectorType >::AdditionalData::absolute_tolerance

Absolute error tolerance for adaptive time stepping.

Definition at line 483 of file ida.h.

◆ relative_tolerance

template<typename VectorType = Vector<double>>
double SUNDIALS::IDA< VectorType >::AdditionalData::relative_tolerance

Relative error tolerance for adaptive time stepping.

Definition at line 488 of file ida.h.

◆ maximum_order

template<typename VectorType = Vector<double>>
unsigned int SUNDIALS::IDA< VectorType >::AdditionalData::maximum_order

Maximum order of BDF.

Definition at line 493 of file ida.h.

◆ output_period

template<typename VectorType = Vector<double>>
double SUNDIALS::IDA< VectorType >::AdditionalData::output_period

Time period between each output.

Definition at line 498 of file ida.h.

◆ ignore_algebraic_terms_for_errors

template<typename VectorType = Vector<double>>
bool SUNDIALS::IDA< VectorType >::AdditionalData::ignore_algebraic_terms_for_errors

Ignore algebraic terms for errors.

Definition at line 503 of file ida.h.

◆ ic_type

template<typename VectorType = Vector<double>>
InitialConditionCorrection SUNDIALS::IDA< VectorType >::AdditionalData::ic_type

Type of correction for initial conditions.

If you do not provide consistent initial conditions, (i.e., conditions for which \(F(y_dot(0), y(0), 0) = 0\)), you can ask SUNDIALS to compute initial conditions for you by using the ic_type parameter at construction time.

Notice that you could in principle use this capabilities to solve for steady state problems by setting y_dot to zero, and asking to compute \(y(0)\) that satisfies \(F(0, y(0), 0) = 0\), however the nonlinear solver used inside IDA may not be robust enough for complex problems with several millions unknowns.

Definition at line 519 of file ida.h.

◆ reset_type

template<typename VectorType = Vector<double>>
InitialConditionCorrection SUNDIALS::IDA< VectorType >::AdditionalData::reset_type

Type of correction for initial conditions to be used after a solver restart.

If you do not have consistent initial conditions after a restart, (i.e., conditions for which F(y_dot(t_restart), y(t_restart), t_restart) = 0), you can ask SUNDIALS to compute the new initial conditions for you by using the reset_type parameter at construction time.

Definition at line 531 of file ida.h.

◆ maximum_non_linear_iterations_ic

template<typename VectorType = Vector<double>>
unsigned SUNDIALS::IDA< VectorType >::AdditionalData::maximum_non_linear_iterations_ic

Maximum number of iterations for Newton method in IC calculation.

Definition at line 536 of file ida.h.

◆ maximum_non_linear_iterations

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

Maximum number of iterations for Newton method during time advancement.

Definition at line 541 of file ida.h.

◆ ls_norm_factor

template<typename VectorType = Vector<double>>
double SUNDIALS::IDA< VectorType >::AdditionalData::ls_norm_factor

Factor to use when converting from the integrator tolerance to the linear solver tolerance

Definition at line 547 of file ida.h.


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