Class SundialsIntegrator#
Defined in File Integrator.h
Nested Relationships#
Nested Types#
Inheritance Relationships#
Base Type#
public Karana::Integrators::Integrator(Class Integrator)
Derived Types#
public Karana::Integrators::ArkExplicitIntegrator(Class ArkExplicitIntegrator)public Karana::Integrators::CVodeIntegrator(Class CVodeIntegrator)
Class Documentation#
-
class SundialsIntegrator : public Karana::Integrators::Integrator#
Abstract base class for Sundials-related integrators.
This class helps manage shared memory, options, and vars. It defines several pure virtual methods for the subclass to override and satisfies all virtual functions of an Integrator.
Subclassed by Karana::Integrators::ArkExplicitIntegrator, Karana::Integrators::CVodeIntegrator
Public Types
Public Functions
-
SundialsIntegrator(const IntegratorType &integrator_type, size_t nstates, DerivativeFunction f, const SundialsIntegrator::IntegratorOptions *options = nullptr)#
Constructor.
- Parameters:
nstates – The number of states to integrate.
f – The derivative function to calculate the derivative of the states.
integrator_type – The type of sundials integrator to use.
options – The integrator options.
-
~SundialsIntegrator()#
-
void setRtol(double rtol)#
Set the relative tolerance used for controlling error vs runtime.
- Parameters:
rtol – The rtol to use for all state variables.
-
void setAtol(double atol)#
Set the scalar absolute tolerance used for controlling error vs runtime.
- Parameters:
atol – The scalar atol to use for all state variables.
-
void setAtol(km::Vec atol_vec)#
Set the vector absolute tolerance used for controlling error vs runtime.
- Parameters:
atol_vec – The vector atol to use; must have same length as state vector.
-
virtual km::Vec getErrorContributions() const = 0#
Get the error contribution for each state.
- Returns:
A vector with the error contribution for each state.
-
virtual void hardReset(size_t nstates) override#
Reset the state size.
- Parameters:
nstates – the new state size
-
virtual void softReset(const km::Ktime &new_t, const km::Vec &new_x) override#
Reset the cached time and state and to the new values going forward.
- Parameters:
new_t – the new time
new_x – the new state
-
virtual km::Ktime advanceTo(const km::Ktime &t_end) override#
Method to advance the system state to a specified time.
- Parameters:
to_t – the desired end time
- Returns:
the actual end time
-
virtual void updateJacobian(const km::Mat &mat) override#
Update the stored Jacobian matrix with the provided one.
This is currently an experimental feature - and not all integrators support the use of Jacobians.
- Parameters:
mat – the new Jacobian matrix
Protected Functions
-
void _initializeCommonMemory()#
Initialize the Sundials integrator.
Called by the constructor, and whenever states are resized. Also calls child class-specific implementations.
-
void _freeCommonMemory()#
Free the memory used by all Sundials integrators.
-
virtual void _initialize() = 0#
Initialize the implementation-specific data structures.
-
virtual void _freeMemory() = 0#
Free all memory, both common and implementation-specific. Must call _freeCommonMemory() to clear shared data structures.
-
virtual void _applyTolerances() = 0#
Copy tolerances from internal options into implementation data structures.
-
virtual void _reInit() = 0#
Reinitialize internal implementation state for a soft reset.
-
virtual void _advanceTo(double t_end) = 0#
Advances implementation to a target time; called by parent class.
- Parameters:
t_end –
Protected Attributes
-
SUNContext _sunctx#
Sundials context object.
-
N_Vector _y#
The current SUNDIALS integrator state.
-
SUNMatrix _a_mat#
Matrix used in the linear solver.
-
SUNLinearSolver _ls#
Linear solver.
-
SUNNonlinearSolver _nls#
Nonlinear solver.
-
km::Vec _dx#
The derivative of the Karana states.
-
km::Mat _jacobian_buffer#
buffer for the Karana Jacobian matrix
Protected Static Functions
-
static int _ode_rhs(sunrealtype t, N_Vector y, N_Vector ydot, void *user_data)#
Calls the derivative function on the provided integrator state, and puts the derivative of the state into ydot.
- Parameters:
t – The time to use when calculating the derivative.
y – The current integrator state.
ydot – The vector that will get populated with the derivative of the state.
user_data – Data used to by the derivative function.
- Returns:
An error code indicating if the call was successful or not.
-
static int _ode_jacobian(sunrealtype t, N_Vector y, N_Vector fy, SUNMatrix J, void *user_data, N_Vector tmp1, N_Vector tmp2, N_Vector tmp3)#
Calculate the jacobian for the nonlinear solver used in ODE solves.
type CVLsJacFn. Invoked by the integrator to get the Jacobian matrix (experimental)
- Parameters:
t – The time to calculate the jacobian at.
y – The state to calculate the jacobian at.
fy – Unused.
J – The matrix that will be populated with the jacobian.
user_data – Any extra data needed to calculate the jacobian.
tmp1 – Unused.
tmp2 – Unused.
tmp3 – Unused.
- Returns:
An error code indicated success or failure of the call.
-
struct IntegratorOptions : public Karana::Integrators::Integrator::IntegratorOptions#
Struct with options for this integrator
Subclassed by Karana::Integrators::ArkExplicitIntegrator::IntegratorOptions, Karana::Integrators::CVodeIntegrator::IntegratorOptions
Public Members
-
double rtol = 1e-4#
The desired relative tolerance used to control the integrator error.
-
double atol = 1e-8#
The desired absolute tolerance used to control the integrator error.
-
km::Vec atol_vec#
The desired absolute tolerance vector used to control the integrator error. This can be used rather than a scalar tolerance to set different values for different states.
-
long max_num_steps = 1e6#
The maximum number of steps the solver can take before throwing a too much work error.
-
double min_step_size = 1e-12#
Minimum acceptable step size before we error out. Defaults to 1ns.
-
double anderson_damping = 1#
Anderson damping in (0, 1] used for implicit integrators with fixed-point solves.
This is used to potentially speed up internal iterations by reusing past iterates (within the same time step). A value of 1 indicates no damping, and smaller values imply a “stronger” damping.
-
int anderson_length = 0#
Number of past iterates to consider with Anderson acceleration.
-
unsigned int max_nl_iters = 3#
Maximum number of nonlinear iterations before we try a smaller step size.
-
NLSOLVER nlsolver = NLSOLVER::FIXEDPOINT#
The Nonlinear solver type to use.
-
bool use_jacobian = false#
set to true to have the integrator use Jacobians (experimental)
-
double rtol = 1e-4#
-
SundialsIntegrator(const IntegratorType &integrator_type, size_t nstates, DerivativeFunction f, const SundialsIntegrator::IntegratorOptions *options = nullptr)#