Class Solver
Defined in File solver.h
Inheritance Relationships
Derived Types
public amici::CVodeSolver(Class CVodeSolver)public amici::IDASolver(Class IDASolver)
Class Documentation
-
class Solver
The Solver class provides a generic interface to CVODES and IDAS solvers, individual realizations are realized in the CVodeSolver and the IDASolver class. All transient private/protected members (CVODES/IDAS memory, interface variables and status flags) are specified as mutable and not included in serialization or equality checks. No solver setting parameter should be marked mutable.
Subclassed by amici::CVodeSolver, amici::IDASolver
Public Types
-
using user_data_type = std::pair<Model*, Solver const*>
Type of what is passed to Sundials solvers as user_data
-
using free_solver_ptr = std::function<void(void*)>
Type of the function to free a raw sundials solver pointer
Public Functions
-
Solver() = default
Default constructor.
-
virtual ~Solver() = default
-
SUNContext get_sun_context() const
Get SUNDIALS context.
- Returns:
context
-
virtual std::string get_class_name() const = 0
Get the name of this class.
- Returns:
Class name.
-
int run(realtype tout) const
runs a forward simulation until the specified timepoint
- Parameters:
tout – next timepoint
- Returns:
status flag
-
int step(realtype tout) const
makes a single step in the simulation
- Parameters:
tout – next timepoint
- Returns:
status flag
-
void run_b(realtype tout) const
runs a backward simulation until the specified timepoint
- Parameters:
tout – next timepoint
-
void setup(realtype t0, Model *model, AmiVector const &x0, AmiVector const &dx0, AmiVectorArray const &sx0, AmiVectorArray const &sdx0) const
Initializes the ami memory object and applies specified options.
- Parameters:
t0 – initial timepoint
model – pointer to the model instance
x0 – initial states
dx0 – initial derivative states
sx0 – initial state sensitivities
sdx0 – initial derivative state sensitivities
-
void setup_b(int *which, realtype tf, Model *model, AmiVector const &xB0, AmiVector const &dxB0, AmiVector const &xQB0) const
Initializes the AMI memory object for the backwards problem.
- Parameters:
which – index of the backward problem, will be set by this routine
tf – final timepoint (initial timepoint for the bwd problem)
model – pointer to the model instance
xB0 – initial adjoint states
dxB0 – initial adjoint derivative states
xQB0 – initial adjoint quadratures
-
void setup_steady_state(realtype t0, Model *model, AmiVector const &x0, AmiVector const &dx0, AmiVector const &xB0, AmiVector const &dxB0, AmiVector const &xQ0) const
Initializes the ami memory for quadrature computation.
- Parameters:
t0 – initial timepoint
model – pointer to the model instance
x0 – initial states
dx0 – initial derivative states
xB0 – initial adjoint states
dxB0 – initial derivative adjoint states
xQ0 – initial quadrature vector
-
void update_and_reinit_states_and_sensitivities(Model *model) const
Reinitializes state and respective sensitivities (if necessary) according to changes in fixedParameters.
- Parameters:
model – pointer to the model instance
-
virtual void get_root_info(int *rootsfound) const = 0
get_root_info extracts information which event occurred
- Parameters:
rootsfound – array with flags indicating whether the respective event occurred
-
virtual void calc_ic(realtype tout1) const = 0
Calculates consistent initial conditions, assumes initial states to be correct (DAE only)
- Parameters:
tout1 – next timepoint to be computed (sets timescale)
-
virtual void calc_ic_b(int which, realtype tout1) const = 0
Calculates consistent initial conditions for the backwards problem, assumes initial states to be correct (DAE only)
- Parameters:
which – identifier of the backwards problem
tout1 – next timepoint to be computed (sets timescale)
-
virtual void solve_b(realtype tBout, int itaskB) const = 0
Solves the backward problem until a predefined timepoint (adjoint only)
- Parameters:
tBout – timepoint until which simulation should be performed
itaskB – task identifier, can be CV_NORMAL or CV_ONE_STEP
-
virtual void turn_off_root_finding() const = 0
Disable rootfinding.
-
SensitivityMethod get_sensitivity_method() const
Return current sensitivity method.
- Returns:
method enum
-
void set_sensitivity_method(SensitivityMethod sensi_meth)
Set sensitivity method.
- Parameters:
sensi_meth –
-
SensitivityMethod get_sensitivity_method_pre_equilibration() const
Return current sensitivity method during preequilibration.
- Returns:
method enum
-
void set_sensitivity_method_pre_equilibration(SensitivityMethod sensi_meth_preeq)
Set sensitivity method for preequilibration.
- Parameters:
sensi_meth_preeq –
-
void switch_forward_sensis_off() const
Disable forward sensitivity integration (used in steady state sim)
-
int get_newton_max_steps() const
Get maximum number of allowed Newton steps for steady state computation.
- Returns:
-
void set_newton_max_steps(int newton_maxsteps)
Set maximum number of allowed Newton steps for steady state computation.
- Parameters:
newton_maxsteps –
-
NewtonDampingFactorMode get_newton_damping_factor_mode() const
Get a state of the damping factor used in the Newton solver.
- Returns:
-
void set_newton_damping_factor_mode(NewtonDampingFactorMode dampingFactorMode)
Turn on/off a damping factor in the Newton method.
- Parameters:
dampingFactorMode –
-
double get_newton_damping_factor_lower_bound() const
Get a lower bound of the damping factor used in the Newton solver.
- Returns:
-
void set_newton_damping_factor_lower_bound(double dampingFactorLowerBound)
Set a lower bound of the damping factor in the Newton solver.
- Parameters:
dampingFactorLowerBound –
-
SensitivityOrder get_sensitivity_order() const
Get sensitivity order.
- Returns:
sensitivity order
-
void set_sensitivity_order(SensitivityOrder sensi)
Set the sensitivity order.
- Parameters:
sensi – sensitivity order
-
double get_relative_tolerance() const
Get the relative tolerances for the forward problem.
Same tolerance is used for the backward problem if not specified differently via setRelativeToleranceASA.
- Returns:
relative tolerances
-
void set_relative_tolerance(double rtol)
Sets the relative tolerances for the forward problem.
Same tolerance is used for the backward problem if not specified differently via setRelativeToleranceASA.
- Parameters:
rtol – relative tolerance (non-negative number)
-
double get_absolute_tolerance() const
Get the absolute tolerances for the forward problem.
Same tolerance is used for the backward problem if not specified differently via setAbsoluteToleranceASA.
- Returns:
absolute tolerances
-
void set_absolute_tolerance(double atol)
Sets the absolute tolerances for the forward problem.
Same tolerance is used for the backward problem if not specified differently via setAbsoluteToleranceASA.
- Parameters:
atol – absolute tolerance (non-negative number)
-
double get_relative_tolerance_fsa() const
Returns the relative tolerances for the forward sensitivity problem.
- Returns:
relative tolerances
-
void set_relative_tolerance_fsa(double rtol)
Sets the relative tolerances for the forward sensitivity problem.
- Parameters:
rtol – relative tolerance (non-negative number)
-
double get_absolute_tolerance_fsa() const
Returns the absolute tolerances for the forward sensitivity problem.
- Returns:
absolute tolerances
-
void set_absolute_tolerance_fsa(double atol)
Sets the absolute tolerances for the forward sensitivity problem.
- Parameters:
atol – absolute tolerance (non-negative number)
-
double get_relative_tolerance_b() const
Returns the relative tolerances for the adjoint sensitivity problem.
- Returns:
relative tolerances
-
void set_relative_tolerance_b(double rtol)
Sets the relative tolerances for the adjoint sensitivity problem.
- Parameters:
rtol – relative tolerance (non-negative number)
-
double get_absolute_tolerance_b() const
Returns the absolute tolerances for the backward problem for adjoint sensitivity analysis.
- Returns:
absolute tolerances
-
void set_absolute_tolerance_b(double atol)
Sets the absolute tolerances for the backward problem for adjoint sensitivity analysis.
- Parameters:
atol – absolute tolerance (non-negative number)
-
double get_relative_tolerance_quadratures() const
Returns the relative tolerance for the quadrature problem.
- Returns:
relative tolerance
-
void set_relative_tolerance_quadratures(double rtol)
sets the relative tolerance for the quadrature problem
- Parameters:
rtol – relative tolerance (non-negative number)
-
double get_absolute_tolerance_quadratures() const
returns the absolute tolerance for the quadrature problem
- Returns:
absolute tolerance
-
void set_absolute_tolerance_quadratures(double atol)
sets the absolute tolerance for the quadrature problem
- Parameters:
atol – absolute tolerance (non-negative number)
-
double get_steady_state_tolerance_factor() const
returns the steady state simulation tolerance factor.
Steady state simulation tolerances are the product of the simulation tolerances and this factor, unless manually set with
set(Absolute/Relative)ToleranceSteadyState().- Returns:
steady state simulation tolerance factor
-
void set_steady_state_tolerance_factor(double factor)
set the steady state simulation tolerance factor.
Steady state simulation tolerances are the product of the simulation tolerances and this factor, unless manually set with
set(Absolute/Relative)ToleranceSteadyState().- Parameters:
factor – tolerance factor (non-negative number)
-
double get_relative_tolerance_steady_state() const
returns the relative tolerance for the steady state problem
- Returns:
relative tolerance
-
void set_relative_tolerance_steady_state(double rtol)
sets the relative tolerance for the steady state problem
- Parameters:
rtol – relative tolerance (non-negative number)
-
double get_absolute_tolerance_steady_state() const
returns the absolute tolerance for the steady state problem
- Returns:
absolute tolerance
-
void set_absolute_tolerance_steady_state(double atol)
sets the absolute tolerance for the steady state problem
- Parameters:
atol – absolute tolerance (non-negative number)
-
double get_steady_state_sensi_tolerance_factor() const
returns the steady state sensitivity simulation tolerance factor.
Steady state sensitivity simulation tolerances are the product of the sensitivity simulation tolerances and this factor, unless manually set with
set(Absolute/Relative)ToleranceSteadyStateSensi().- Returns:
steady state simulation tolerance factor
-
void set_steady_state_sensi_tolerance_factor(double factor)
set the steady state sensitivity simulation tolerance factor.
Steady state sensitivity simulation tolerances are the product of the sensitivity simulation tolerances and this factor, unless manually set with
set(Absolute/Relative)ToleranceSteadyStateSensi().- Parameters:
factor – tolerance factor (non-negative number)
-
double get_relative_tolerance_steady_state_sensi() const
returns the relative tolerance for the sensitivities of the steady state problem
- Returns:
relative tolerance
-
void set_relative_tolerance_steady_state_sensi(double rtol)
sets the relative tolerance for the sensitivities of the steady state problem
- Parameters:
rtol – relative tolerance (non-negative number)
-
double get_absolute_tolerance_steady_state_sensi() const
returns the absolute tolerance for the sensitivities of the steady state problem
- Returns:
absolute tolerance
-
void set_absolute_tolerance_steady_state_sensi(double atol)
sets the absolute tolerance for the sensitivities of the steady state problem
- Parameters:
atol – absolute tolerance (non-negative number)
-
long int get_max_steps() const
returns the maximum number of solver steps for the forward problem
- Returns:
maximum number of solver steps
-
void set_max_steps(long int maxsteps)
sets the maximum number of solver steps for the forward problem
- Parameters:
maxsteps – maximum number of solver steps (positive number)
-
double get_max_time() const
Returns the maximum time allowed for integration.
- Returns:
Time in seconds
-
void set_max_time(double maxtime)
Set the maximum CPU time allowed for integration.
- Parameters:
maxtime – Time in seconds. Zero means infinite time.
-
void start_timer() const
Start timer for tracking integration time.
-
bool time_exceeded(int interval = 1) const
Check whether maximum integration time was exceeded.
- Parameters:
interval – Only check the time every
intervalths call to avoid potentially relatively expensive syscalls- Returns:
True if the maximum integration time was exceeded, false otherwise.
-
long int get_max_steps_backward_problem() const
returns the maximum number of solver steps for the backward problem
- Returns:
maximum number of solver steps
-
void set_max_steps_backward_problem(long int maxsteps)
sets the maximum number of solver steps for the backward problem
Note
default behaviour (100 times the value for the forward problem) can be restored by passing maxsteps=0
- Parameters:
maxsteps – maximum number of solver steps (non-negative number)
-
LinearMultistepMethod get_linear_multistep_method() const
returns the linear system multistep method
- Returns:
linear system multistep method
-
void set_linear_multistep_method(LinearMultistepMethod lmm)
sets the linear system multistep method
- Parameters:
lmm – linear system multistep method
-
NonlinearSolverIteration get_non_linear_solver_iteration() const
returns the nonlinear system solution method
- Returns:
-
void set_non_linear_solver_iteration(NonlinearSolverIteration iter)
sets the nonlinear system solution method
- Parameters:
iter – nonlinear system solution method
-
InterpolationType get_interpolation_type() const
get_interpolation_type
- Returns:
-
void set_interpolation_type(InterpolationType interpType)
sets the interpolation of the forward solution that is used for the backwards problem
- Parameters:
interpType – interpolation type
-
int get_state_ordering() const
Gets KLU / SuperLUMT state ordering mode.
- Returns:
State-ordering as integer according to SUNLinSolKLU::StateOrdering or SUNLinSolSuperLUMT::StateOrdering (which differ).
-
void set_state_ordering(int ordering)
Sets KLU / SuperLUMT state ordering mode.
This only applies when linsol is set to LinearSolver::KLU or LinearSolver::SuperLUMT. Mind the difference between SUNLinSolKLU::StateOrdering and SUNLinSolSuperLUMT::StateOrdering.
- Parameters:
ordering – state ordering
-
bool get_stability_limit_flag() const
returns stability limit detection mode
- Returns:
stldet can be false (deactivated) or true (activated)
-
void set_stability_limit_flag(bool stldet)
set stability limit detection mode
- Parameters:
stldet – can be false (deactivated) or true (activated)
-
LinearSolver get_linear_solver() const
get_linear_solver
- Returns:
-
void set_linear_solver(LinearSolver linsol)
set_linear_solver
- Parameters:
linsol –
-
InternalSensitivityMethod get_internal_sensitivity_method() const
returns the internal sensitivity method
- Returns:
internal sensitivity method
-
void set_internal_sensitivity_method(InternalSensitivityMethod ism)
sets the internal sensitivity method
- Parameters:
ism – internal sensitivity method
-
RDataReporting get_return_data_reporting_mode() const
returns the ReturnData reporting mode
- Returns:
ReturnData reporting mode
-
void set_return_data_reporting_mode(RDataReporting rdrm)
sets the ReturnData reporting mode
- Parameters:
rdrm – ReturnData reporting mode
-
void write_solution(realtype &t, AmiVector &x, AmiVector &dx, AmiVectorArray &sx, AmiVector &xQ) const
write solution from forward simulation
- Parameters:
t – time
x – state
dx – derivative state
sx – state sensitivity
xQ – quadrature
-
void write_solution(realtype &t, AmiVector &x, AmiVector &dx, AmiVectorArray &sx) const
write solution from forward simulation
- Parameters:
t – time
x – state
dx – derivative state
sx – state sensitivity
-
void write_solution(SolutionState &sol) const
write solution from forward simulation
- Parameters:
sol – solution state
-
void write_solution(realtype t, SolutionState &sol) const
write solution from forward simulation
- Parameters:
t – Time for which to retrieve the solution (interpolated if necessary). Must be greater than or equal to the initial timepoint and less than or equal to the current timepoint.
sol – solution state
-
void write_solution_b(realtype &t, AmiVector &xB, AmiVector &dxB, AmiVector &xQB, int which) const
write solution from backward simulation
- Parameters:
t – time
xB – adjoint state
dxB – adjoint derivative state
xQB – adjoint quadrature
which – index of adjoint problem
-
AmiVector const &get_state(realtype t) const
Access state solution at time t.
- Parameters:
t – time
- Returns:
x or interpolated solution dky
-
AmiVector const &get_derivative_state(realtype t) const
Access derivative state solution at time t.
- Parameters:
t – time
- Returns:
dx or interpolated solution dky
-
AmiVectorArray const &get_state_sensitivity(realtype t) const
Access state sensitivity solution at time t.
- Parameters:
t – time
- Returns:
(interpolated) solution sx
-
AmiVector const &get_adjoint_state(int which, realtype t) const
Access adjoint solution at time t.
- Parameters:
which – adjoint problem index
t – time
- Returns:
(interpolated) solution xB
-
AmiVector const &get_adjoint_derivative_state(int which, realtype t) const
Access adjoint derivative solution at time t.
- Parameters:
which – adjoint problem index
t – time
- Returns:
(interpolated) solution dxB
-
AmiVector const &get_adjoint_quadrature(int which, realtype t) const
Access adjoint quadrature solution at time t.
- Parameters:
which – adjoint problem index
t – time
- Returns:
(interpolated) solution xQB
-
AmiVector const &get_quadrature(realtype t) const
Access quadrature solution at time t.
- Parameters:
t – time
- Returns:
(interpolated) solution xQ
-
virtual void reinit(realtype t0, AmiVector const &yy0, AmiVector const &yp0) const = 0
Reinitializes the states in the solver after an event occurrence.
- Parameters:
t0 – reinitialization timepoint
yy0 – initial state variables
yp0 – initial derivative state variables (DAE only)
-
virtual void sens_reinit(AmiVectorArray const &yyS0, AmiVectorArray const &ypS0) const = 0
Reinitializes the state sensitivities in the solver after an event occurrence.
- Parameters:
yyS0 – new state sensitivity
ypS0 – new derivative state sensitivities (DAE only)
-
virtual void sens_toggle_off() const = 0
Switches off computation of state sensitivities without deallocating the memory for sensitivities.
-
virtual void reinit_b(int which, realtype tB0, AmiVector const &yyB0, AmiVector const &ypB0) const = 0
Reinitializes the adjoint states after an event occurrence.
- Parameters:
which – identifier of the backwards problem
tB0 – reinitialization timepoint
yyB0 – new adjoint state
ypB0 – new adjoint derivative state
-
virtual void reinit_quad_b(int which, AmiVector const &yQB0) const = 0
Reinitialize the adjoint states after an event occurrence.
- Parameters:
which – identifier of the backwards problem
yQB0 – new adjoint quadrature state
-
realtype get_cpu_time_b() const
Reads out the CPU time needed for backward solve.
- Returns:
cpu_timeB
-
int nx() const
number of states with which the solver was initialized
- Returns:
x.getLength()
-
int nplist() const
number of parameters with which the solver was initialized
- Returns:
sx.getLength()
-
int nquad() const
number of quadratures with which the solver was initialized
- Returns:
xQB.getLength()
-
inline bool computing_fsa() const
check if FSA is being computed
- Returns:
flag
-
inline bool computing_asa() const
check if ASA is being computed
- Returns:
flag
-
void reset_diagnosis() const
Resets vectors containing diagnosis information.
-
void store_diagnosis() const
Stores diagnosis information from solver memory block for forward problem.
-
void store_diagnosis_b(int which) const
Stores diagnosis information from solver memory block for backward problem.
- Parameters:
which – identifier of the backwards problem
-
inline std::vector<int> const &get_num_steps() const
Accessor ns.
- Returns:
ns
-
inline std::vector<int> const &get_num_steps_b() const
Accessor nsB.
- Returns:
nsB
-
inline std::vector<int> const &get_num_rhs_evals() const
Accessor nrhs.
- Returns:
nrhs
-
inline std::vector<int> const &get_num_rhs_evals_b() const
Accessor nrhsB.
- Returns:
nrhsB
-
inline std::vector<int> const &get_num_err_test_fails() const
Accessor netf.
- Returns:
netf
-
inline std::vector<int> const &get_num_err_test_fails_b() const
Accessor netfB.
- Returns:
netfB
-
inline std::vector<int> const &get_num_non_lin_solv_conv_fails() const
Accessor nnlscf.
- Returns:
nnlscf
-
inline std::vector<int> const &get_num_non_lin_solv_conv_fails_b() const
Accessor nnlscfB.
- Returns:
nnlscfB
-
inline std::vector<int> const &get_last_order() const
Accessor order.
- Returns:
order
-
inline bool get_newton_step_steady_state_check() const
Returns how convergence checks for steadystate computation are performed. If activated, convergence checks are limited to every 25 steps in the simulation solver to limit performance impact.
- Returns:
boolean flag indicating newton step (true) or the right hand side (false)
-
inline bool get_sensi_steady_state_check() const
Returns how convergence checks for steadystate computation are performed.
- Returns:
boolean flag indicating state and sensitivity equations (true) or only state variables (false).
-
inline void set_newton_step_steady_state_check(bool const flag)
Sets how convergence checks for steadystate computation are performed.
- Parameters:
flag – boolean flag to pick newton step (true) or the right hand side (false, default)
-
inline void set_sensi_steady_state_check(bool const flag)
Sets for which variables convergence checks for steadystate computation are performed.
- Parameters:
flag – boolean flag to pick state and sensitivity equations (true, default) or only state variables (false).
-
void set_max_nonlin_iters(int max_nonlin_iters)
Set the maximum number of nonlinear solver iterations permitted per step.
- Parameters:
max_nonlin_iters – maximum number of nonlinear solver iterations
-
int get_max_nonlin_iters() const
Get the maximum number of nonlinear solver iterations permitted per step.
- Returns:
maximum number of nonlinear solver iterations
-
void set_max_conv_fails(int max_conv_fails)
Set the maximum number of nonlinear solver convergence failures permitted per step.
- Parameters:
max_conv_fails – maximum number of nonlinear solver convergence
-
int get_max_conv_fails() const
Get the maximum number of nonlinear solver convergence failures permitted per step.
- Returns:
maximum number of nonlinear solver convergence
-
void set_constraints(std::vector<realtype> const &constraints)
Set constraints on the model state.
See https://sundials.readthedocs.io/en/latest/cvode/Usage/index.html#c.CVodeSetConstraints.
- Parameters:
constraints –
-
inline std::vector<realtype> get_constraints() const
Get constraints on the model state.
- Returns:
constraints
Protected Functions
-
virtual void set_stop_time(realtype tstop) const = 0
Sets a timepoint at which the simulation will be stopped.
- Parameters:
tstop – timepoint until which simulation should be performed
-
virtual int solve(realtype tout, int itask) const = 0
Solves the forward problem until a predefined timepoint.
- Parameters:
tout – timepoint until which simulation should be performed
itask – task identifier, can be CV_NORMAL or CV_ONE_STEP
- Returns:
status flag indicating success of execution
-
virtual int solve_f(realtype tout, int itask, int *ncheckPtr) const = 0
Solves the forward problem until a predefined timepoint (adjoint only)
- Parameters:
tout – timepoint until which simulation should be performed
itask – task identifier, can be CV_NORMAL or CV_ONE_STEP
ncheckPtr – pointer to a number that counts the internal checkpoints
- Returns:
status flag indicating success of execution
-
virtual void reinit_post_process_f(realtype tnext) const = 0
reinit_post_process_f postprocessing of the solver memory after a discontinuity in the forward problem
- Parameters:
tnext – next timepoint (defines integration direction)
-
virtual void reinit_post_process_b(realtype tnext) const = 0
reinit_post_process_b postprocessing of the solver memory after a discontinuity in the backward problem
- Parameters:
tnext – next timepoint (defines integration direction)
-
virtual void get_sens() const = 0
extracts the state sensitivity at the current timepoint from solver memory and writes it to the sx member variable
-
virtual void get_b(int which) const = 0
extracts the adjoint state at the current timepoint from solver memory and writes it to the xB member variable
- Parameters:
which – index of the backwards problem
-
virtual void get_quad_b(int which) const = 0
extracts the adjoint quadrature state at the current timepoint from solver memory and writes it to the xQB member variable
- Parameters:
which – index of the backwards problem
-
virtual void get_quad(realtype &t) const = 0
extracts the quadrature at the current timepoint from solver memory and writes it to the xQ member variable
- Parameters:
t – timepoint for quadrature extraction
-
virtual void init(realtype t0, AmiVector const &x0, AmiVector const &dx0) const = 0
Initializes the states at the specified initial timepoint.
- Parameters:
t0 – initial timepoint
x0 – initial states
dx0 – initial derivative states
-
virtual void init_steady_state(realtype t0, AmiVector const &x0, AmiVector const &dx0) const = 0
Initializes the states at the specified initial timepoint.
- Parameters:
t0 – initial timepoint
x0 – initial states
dx0 – initial derivative states
-
virtual void sens_init_1(AmiVectorArray const &sx0, AmiVectorArray const &sdx0) const = 0
Initializes the forward sensitivities.
- Parameters:
sx0 – initial states sensitivities
sdx0 – initial derivative states sensitivities
-
virtual void b_init(int which, realtype tf, AmiVector const &xB0, AmiVector const &dxB0) const = 0
Initialize the adjoint states at the specified final timepoint.
- Parameters:
which – identifier of the backwards problem
tf – final timepoint
xB0 – initial adjoint state
dxB0 – initial adjoint derivative state
-
virtual void qb_init(int which, AmiVector const &xQB0) const = 0
Initialize the quadrature states at the specified final timepoint.
- Parameters:
which – identifier of the backwards problem
xQB0 – initial adjoint quadrature state
-
virtual void root_init(int ne) const = 0
Initializes the root-finding for events.
- Parameters:
ne – number of different events
-
void initialize_non_linear_solver_sens(Model const *model) const
Initialize non-linear solver for sensitivities.
- Parameters:
model – Model instance
-
virtual void set_dense_jac_fn() const = 0
Set the dense Jacobian function.
-
virtual void set_sparse_jac_fn() const = 0
sets the sparse Jacobian function
-
virtual void set_band_jac_fn() const = 0
sets the banded Jacobian function
-
virtual void set_jac_times_vec_fn() const = 0
sets the Jacobian vector multiplication function
-
virtual void set_dense_jac_fn_b(int which) const = 0
sets the dense Jacobian function
- Parameters:
which – identifier of the backwards problem
-
virtual void set_sparse_jac_fn_b(int which) const = 0
sets the sparse Jacobian function
- Parameters:
which – identifier of the backwards problem
-
virtual void set_band_jac_fn_b(int which) const = 0
sets the banded Jacobian function
- Parameters:
which – identifier of the backwards problem
-
virtual void set_jac_times_vec_fn_b(int which) const = 0
sets the Jacobian vector multiplication function
- Parameters:
which – identifier of the backwards problem
-
virtual void set_sparse_jac_fn_ss() const = 0
sets the sparse Jacobian function for backward steady state case
-
virtual void allocate_solver() const = 0
Create specifies solver method and initializes solver memory for the forward problem.
-
virtual void set_ss_tolerances(double rtol, double atol) const = 0
sets scalar relative and absolute tolerances for the forward problem
- Parameters:
rtol – relative tolerances
atol – absolute tolerances
-
virtual void set_sens_ss_tolerances(double rtol, double const *atol) const = 0
activates sets scalar relative and absolute tolerances for the sensitivity variables
- Parameters:
rtol – relative tolerances
atol – array of absolute tolerances for every sensitivity variable
-
virtual void set_sens_err_con(bool error_corr) const = 0
SetSensErrCon specifies whether error control is also enforced for sensitivities for the forward problem
- Parameters:
error_corr – activation flag
-
virtual void set_quad_err_con_b(int which, bool flag) const = 0
Specifies whether error control is also enforced for the backward quadrature problem.
- Parameters:
which – identifier of the backwards problem
flag – activation flag
-
virtual void set_quad_err_con(bool flag) const = 0
Specifies whether error control is also enforced for the forward quadrature problem.
- Parameters:
flag – activation flag
-
virtual void set_err_handler_fn() const
Attaches the error handler function to the solver.
-
virtual void set_user_data() const = 0
Attaches the user data to the forward problem.
-
virtual void set_user_data_b(int which) const = 0
attaches the user data to the backward problem
- Parameters:
which – identifier of the backwards problem
-
virtual void set_max_num_steps(long int mxsteps) const = 0
specifies the maximum number of steps for the forward problem
Note
in contrast to the SUNDIALS method, this sets the overall maximum, not the maximum between output times.
- Parameters:
mxsteps – number of steps
-
virtual void set_max_num_steps_b(int which, long int mxstepsB) const = 0
specifies the maximum number of steps for the forward problem
Note
in contrast to the SUNDIALS method, this sets the overall maximum, not the maximum between output times.
- Parameters:
which – identifier of the backwards problem
mxstepsB – number of steps
-
virtual void set_stab_lim_det(int stldet) const = 0
activates stability limit detection for the forward problem
- Parameters:
stldet – flag for stability limit detection (TRUE or FALSE)
-
virtual void set_stab_lim_det_b(int which, int stldet) const = 0
activates stability limit detection for the backward problem
- Parameters:
which – identifier of the backwards problem
stldet – flag for stability limit detection (TRUE or FALSE)
-
virtual void set_id(Model const *model) const = 0
specify algebraic/differential components (DAE only)
- Parameters:
model – model specification
-
virtual void set_suppress_alg(bool flag) const = 0
deactivates error control for algebraic components (DAE only)
- Parameters:
flag – deactivation flag
-
virtual void set_sens_params(realtype const *p, realtype const *pbar, int const *plist) const = 0
specifies the scaling and indexes for sensitivity computation
- Parameters:
p – parameters
pbar – parameter scaling constants
plist – parameter index list
-
virtual void get_dky(realtype t, int k) const = 0
interpolates the (derivative of the) solution at the requested timepoint
- Parameters:
t – timepoint
k – derivative order
-
virtual void get_dky_b(realtype t, int k, int which) const = 0
interpolates the (derivative of the) solution at the requested timepoint
- Parameters:
t – timepoint
k – derivative order
which – index of backward problem
-
virtual void get_sens_dky(realtype t, int k) const = 0
interpolates the (derivative of the) solution at the requested timepoint
- Parameters:
t – timepoint
k – derivative order
-
virtual void get_quad_dky_b(realtype t, int k, int which) const = 0
interpolates the (derivative of the) solution at the requested timepoint
- Parameters:
t – timepoint
k – derivative order
which – index of backward problem
-
virtual void get_quad_dky(realtype t, int k) const = 0
interpolates the (derivative of the) solution at the requested timepoint
- Parameters:
t – timepoint
k – derivative order
-
virtual void adj_init() const = 0
initializes the adjoint problem
-
virtual void quad_init(AmiVector const &xQ0) const = 0
initializes the quadratures
- Parameters:
xQ0 – vector with initial values for xQ
-
virtual void allocate_solver_b(int *which) const = 0
Specifies solver method and initializes solver memory for the backward problem.
- Parameters:
which – identifier of the backwards problem
-
virtual void set_ss_tolerances_b(int which, realtype relTolB, realtype absTolB) const = 0
sets relative and absolute tolerances for the backward problem
- Parameters:
which – identifier of the backwards problem
relTolB – relative tolerances
absTolB – absolute tolerances
-
virtual void quad_ss_tolerances_b(int which, realtype reltolQB, realtype abstolQB) const = 0
sets relative and absolute tolerances for the quadrature backward problem
- Parameters:
which – identifier of the backwards problem
reltolQB – relative tolerances
abstolQB – absolute tolerances
-
virtual void quad_ss_tolerances(realtype reltolQB, realtype abstolQB) const = 0
sets relative and absolute tolerances for the quadrature problem
- Parameters:
reltolQB – relative tolerances
abstolQB – absolute tolerances
-
virtual void get_num_steps(void const *ami_mem, long int *numsteps) const = 0
reports the number of solver steps
- Parameters:
ami_mem – pointer to the solver memory instance (can be from forward or backward problem)
numsteps – output array
-
virtual void get_num_rhs_evals(void const *ami_mem, long int *numrhsevals) const = 0
reports the number of right hand evaluations
- Parameters:
ami_mem – pointer to the solver memory instance (can be from forward or backward problem)
numrhsevals – output array
-
virtual void get_num_err_test_fails(void const *ami_mem, long int *numerrtestfails) const = 0
reports the number of local error test failures
- Parameters:
ami_mem – pointer to the solver memory instance (can be from forward or backward problem)
numerrtestfails – output array
-
virtual void get_num_non_lin_solv_conv_fails(void const *ami_mem, long int *numnonlinsolvconvfails) const = 0
reports the number of nonlinear convergence failures
- Parameters:
ami_mem – pointer to the solver memory instance (can be from forward or backward problem)
numnonlinsolvconvfails – output array
-
virtual void get_last_order(void const *ami_mem, int *order) const = 0
Reports the order of the integration method during the last internal step.
- Parameters:
ami_mem – pointer to the solver memory instance (can be from forward or backward problem)
order – output array
-
void initialize_linear_solver(Model const *model) const
Initializes and sets the linear solver for the forward problem.
- Parameters:
model – pointer to the model object
-
void initialize_non_linear_solver() const
Sets the non-linear solver.
-
virtual void set_linear_solver() const = 0
Sets the linear solver for the forward problem.
-
virtual void set_linear_solver_b(int which) const = 0
Sets the linear solver for the backward problem.
- Parameters:
which – index of the backward problem
-
virtual void set_non_linear_solver() const = 0
Set the non-linear solver for the forward problem.
-
virtual void set_non_linear_solver_b(int which) const = 0
Set the non-linear solver for the backward problem.
- Parameters:
which – index of the backward problem
-
virtual void set_non_linear_solver_sens() const = 0
Set the non-linear solver for sensitivities.
-
void initialize_linear_solver_b(Model const *model, int which) const
Initializes the linear solver for the backward problem.
- Parameters:
model – pointer to the model object
which – index of the backward problem
-
void initialize_non_linear_solver_b(int which) const
Initializes the non-linear solver for the backward problem.
- Parameters:
which – index of the backward problem
-
virtual Model const *get_model() const = 0
Accessor function to the model stored in the user data
- Returns:
user data model
-
bool get_init_done() const
checks whether memory for the forward problem has been allocated
- Returns:
proxy for solverMemory->(cv|ida)_MallocDone
-
bool get_sens_init_done() const
checks whether memory for forward sensitivities has been allocated
- Returns:
proxy for solverMemory->(cv|ida)_SensMallocDone
-
bool get_adj_init_done() const
checks whether memory for forward interpolation has been allocated
- Returns:
proxy for solverMemory->(cv|ida)_adjMallocDone
-
bool get_init_done_b(int which) const
checks whether memory for the backward problem has been allocated
- Parameters:
which – adjoint problem index
- Returns:
proxy for solverMemoryB->(cv|ida)_MallocDone
-
bool get_quad_init_done_b(int which) const
checks whether memory for backward quadratures has been allocated
- Parameters:
which – adjoint problem index
- Returns:
proxy for solverMemoryB->(cv|ida)_QuadMallocDone
-
bool get_quad_init_done() const
checks whether memory for quadratures has been allocated
- Returns:
proxy for solverMemory->(cv|ida)_QuadMallocDone
-
virtual void diag() const = 0
attaches a diagonal linear solver to the forward problem
-
virtual void diag_b(int which) const = 0
attaches a diagonal linear solver to the backward problem
- Parameters:
which – identifier of the backwards problem
-
void reset_mutable_memory(int nx, int nplist, int nquad) const
resets solverMemory and solverMemoryB
- Parameters:
nx – new number of state variables
nplist – new number of sensitivity parameters
nquad – new number of quadratures (only differs from nplist for higher order sensitivity computation)
-
virtual void *get_adj_b_mem(void *ami_mem, int which) const = 0
Retrieves the solver memory instance for the backward problem.
- Parameters:
which – identifier of the backwards problem
ami_mem – pointer to the forward solver memory instance
- Returns:
A (void *) pointer to the CVODES memory allocated for the backward problem.
-
void apply_tolerances() const
updates solver tolerances according to the currently specified member variables
-
void apply_tolerances_fsa() const
updates FSA solver tolerances according to the currently specified member variables
-
void apply_tolerances_asa(int which) const
updates ASA solver tolerances according to the currently specified member variables
- Parameters:
which – identifier of the backwards problem
-
void apply_quad_tolerances_asa(int which) const
updates ASA quadrature solver tolerances according to the currently specified member variables
- Parameters:
which – identifier of the backwards problem
-
void apply_quad_tolerances() const
updates quadrature solver tolerances according to the currently specified member variables
-
void apply_sensitivity_tolerances() const
updates all sensitivity solver tolerances according to the currently specified member variables
-
virtual void apply_constraints() const
Apply the constraints to the solver.
-
void set_init_done() const
sets that memory for the forward problem has been allocated
-
void set_sens_init_done() const
sets that memory for forward sensitivities has been allocated
-
void set_adj_init_done() const
sets that memory for forward interpolation has been allocated
-
void set_init_done_b(int which) const
sets that memory for the backward problem has been allocated
- Parameters:
which – adjoint problem index
-
void set_quad_init_done_b(int which) const
sets that memory for backward quadratures has been allocated
- Parameters:
which – adjoint problem index
-
void set_quad_init_done() const
sets that memory for quadratures has been allocated
-
void check_sensitivity_method(SensitivityMethod sensi_meth, bool preequilibration) const
Sets sensitivity method (for simulation or preequilibration)
- Parameters:
sensi_meth – new value for sensi_meth[_preeq]
preequilibration – flag indicating preequilibration or simulation
-
virtual void apply_max_nonlin_iters() const = 0
Apply the maximum number of nonlinear solver iterations permitted per step.
-
virtual void apply_max_conv_fails() const = 0
Apply the maximum number of nonlinear solver convergence failures permitted per step.
-
virtual void apply_max_step_size() const = 0
Apply the allowed maximum stepsize to the solver.
Protected Attributes
-
sundials::Context sunctx_
SUNDIALS context
-
mutable std::unique_ptr<void, free_solver_ptr> solver_memory_
pointer to solver memory block
-
mutable std::vector<std::unique_ptr<void, free_solver_ptr>> solver_memory_B_
pointer to solver memory block
-
mutable user_data_type user_data_
Sundials user_data_
-
InternalSensitivityMethod ism_ = {InternalSensitivityMethod::simultaneous}
internal sensitivity method flag used to select the sensitivity solution method. Only applies for Forward Sensitivities.
-
LinearMultistepMethod lmm_ = {LinearMultistepMethod::BDF}
specifies the linear multistep method.
-
NonlinearSolverIteration iter_ = {NonlinearSolverIteration::newton}
specifies the type of nonlinear solver iteration
-
InterpolationType interp_type_ = {InterpolationType::polynomial}
interpolation type for the forward problem solution which is then used for the backwards problem.
-
long int maxsteps_ = {10000}
maximum number of allowed integration steps
-
std::chrono::duration<double, std::ratio<1>> maxtime_ = {0}
Maximum CPU-time for integration in seconds
-
mutable std::unique_ptr<SUNLinSolWrapper> linear_solver_
linear solver for the forward problem
-
mutable std::unique_ptr<SUNLinSolWrapper> linear_solver_B_
linear solver for the backward problem
-
mutable std::unique_ptr<SUNNonLinSolWrapper> non_linear_solver_
non-linear solver for the forward problem
-
mutable std::unique_ptr<SUNNonLinSolWrapper> non_linear_solver_B_
non-linear solver for the backward problem
-
mutable std::unique_ptr<SUNNonLinSolWrapper> non_linear_solver_sens_
non-linear solver for the sensitivities
-
mutable bool solver_was_called_F_ = {false}
flag indicating whether the forward solver has been called
-
mutable bool solver_was_called_B_ = {false}
flag indicating whether the backward solver has been called
-
mutable AmiVectorArray sx_ = {0, 0, sunctx_}
state sensitivities interface variable (dimension: nx_solver x nplist)
-
mutable AmiVectorArray sdx_ = {0, 0, sunctx_}
state derivative sensitivities dummy (dimension: nx_solver x nplist)
-
mutable AmiVector xQB_ = {0, sunctx_}
adjoint quadrature interface variable (dimension: nJ x nplist)
-
mutable bool force_reinit_postprocess_F_ = {false}
flag to force reInitPostProcessF before next call to solve
-
mutable bool force_reinit_postprocess_B_ = {false}
flag to force reInitPostProcessB before next call to solveB
-
mutable bool sens_initialized_ = {false}
flag indicating whether sensInit1 was called
-
using user_data_type = std::pair<Model*, Solver const*>