Class CVodeIntegrator#
Defined in File Integrator.h
Nested Relationships#
Nested Types#
Inheritance Relationships#
Base Type#
public Karana::Math::Integrator(Class Integrator)
Class Documentation#
-
class CVodeIntegrator : public Karana::Math::Integrator#
The multistep CVode integrator from the SUNDIALS suite.
Public Types
Public Functions
-
CVodeIntegrator(size_t nstates, DerivativeFunction f, const IntegratorType &integrator_type, const 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 cvode integrator to use.
options – The integrator options.
-
~CVodeIntegrator()#
-
virtual std::string_view typeString() const override#
Return the integrator type as a string.
- Returns:
the integrator type
-
virtual void hardReset(size_t nstates) override#
Reset the state size.
- Parameters:
nstates – the new state size
-
virtual void softReset(const Ktime &new_t, const 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 Ktime advanceTo(const Ktime &t_end) override#
Method to advance the sytem state to a specified time.
- Parameters:
to_t – the desired end time
- Returns:
the actual end time
-
virtual void updateJacobian(const 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
-
virtual std::unique_ptr<Integrator::IntegratorOptions> getOptions() override#
Get the options associated with the integrator.
- Returns:
The options associated with the integrator.
-
virtual IntegratorType getIntegratorType() const override#
Get the IntegratorType for this integrator.
- Returns:
The IntegratorType for this integrator.
-
Vec getErrorContributions() const#
Get the error contributions of the state components.
This is helpful for identifying problematic state components
- Returns:
The error contributions.
-
int getCurrentOrder() const#
Return the current order of the integrator.
- Returns:
The integrator order
-
double getLastStepSize() const#
Return the last step size used.
- Returns:
The last step size
Protected Attributes
-
SUNContext _sunctx = nullptr#
Sundials context object.
-
N_Vector _y = nullptr#
The current integrator state.
-
SUNMatrix _a_mat = nullptr#
Matrix used in the linear solver.
-
SUNLinearSolver _ls = nullptr#
Linear solver.
-
SUNNonlinearSolver _nls = nullptr#
nonLinear solver
-
void *_cvode_mem#
Memory used by the cvode integrator.
-
IntegratorType _integrator_type#
The integrator type.
Protected Static Functions
-
static int _cvode_f(realtype t, N_Vector y, N_Vector ydot, void *user_data)#
Calles 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 deriative function.
- Returns:
An error code indicating if the call was successful or not.
-
static int _cvode_jacobian(realtype 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.
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 willl 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::Math::Integrator::IntegratorOptions#
Struct with options for this integrator
Public Functions
-
IntegratorOptions() = default#
-
IntegratorOptions(const IntegratorOptions &p)#
Copy constructor.
- Parameters:
p – The IntegratorOptions to copy values from.
-
IntegratorOptions &operator=(const IntegratorOptions &p)#
Assignment operator.
- Parameters:
p – The IntegratorOptions to copy values from.
- Returns:
A reference to the newly assigned IntegratorOptions.
-
virtual ~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.
-
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.
-
unsigned int anderson_accel = 0#
Anderson acceleration used with the FixedPoint solver.
-
NLSOLVER nlsolver = NLSOLVER::FIXEDPOINT#
The Nonlinear solver type to use.
-
bool use_jacobian = false#
set to true to have the integrator use Jacobians (experimental)
-
IntegratorOptions() = default#
-
CVodeIntegrator(size_t nstates, DerivativeFunction f, const IntegratorType &integrator_type, const IntegratorOptions *options = nullptr)#