Functions | |
CASADI_EXPORT Function | casadi::integrator (const std::string &name, const std::string &solver, const SXDict &dae, const Dict &opts=Dict()) |
CASADI_EXPORT Function | casadi::integrator (const std::string &name, const std::string &solver, const MXDict &dae, const Dict &opts=Dict()) |
CASADI_EXPORT Function | casadi::integrator (const std::string &name, const std::string &solver, const Function &dae, const Dict &opts=Dict()) |
CASADI_EXPORT Function | casadi::integrator (const std::string &name, const std::string &solver, const SXDict &dae, double t0, const std::vector< double > &tout, const Dict &opts=Dict()) |
CASADI_EXPORT Function | casadi::integrator (const std::string &name, const std::string &solver, const MXDict &dae, double t0, const std::vector< double > &tout, const Dict &opts=Dict()) |
CASADI_EXPORT Function | casadi::integrator (const std::string &name, const std::string &solver, const Function &dae, double t0, const std::vector< double > &tout, const Dict &opts=Dict()) |
CASADI_EXPORT Function | casadi::integrator (const std::string &name, const std::string &solver, const SXDict &dae, double t0, double tf, const Dict &opts=Dict()) |
CASADI_EXPORT Function | casadi::integrator (const std::string &name, const std::string &solver, const MXDict &dae, double t0, double tf, const Dict &opts=Dict()) |
CASADI_EXPORT Function | casadi::integrator (const std::string &name, const std::string &solver, const Function &dae, double t0, double tf, const Dict &opts=Dict()) |
CASADI_EXPORT bool | casadi::has_integrator (const std::string &name) |
Check if a particular plugin is available. More... | |
CASADI_EXPORT void | casadi::load_integrator (const std::string &name) |
Explicitly load a plugin dynamically. More... | |
CASADI_EXPORT std::string | casadi::doc_integrator (const std::string &name) |
Get the documentation string for a plugin. More... | |
CASADI_EXPORT std::vector< std::string > | casadi::integrator_in () |
Get input scheme of integrators. More... | |
CASADI_EXPORT std::vector< std::string > | casadi::integrator_out () |
Get integrator output scheme of integrators. More... | |
CASADI_EXPORT std::string | casadi::integrator_in (casadi_int ind) |
Get integrator input scheme name by index. More... | |
CASADI_EXPORT std::string | casadi::integrator_out (casadi_int ind) |
Get output scheme name by index. More... | |
CASADI_EXPORT casadi_int | casadi::integrator_n_in () |
Get the number of integrator inputs. More... | |
CASADI_EXPORT casadi_int | casadi::integrator_n_out () |
Get the number of integrator outputs. More... | |
CASADI_EXPORT std::vector< std::string > | casadi::dyn_in () |
Get input scheme of simulators. More... | |
CASADI_EXPORT std::vector< std::string > | casadi::dyn_out () |
Get simulator output scheme of simulators. More... | |
CASADI_EXPORT std::string | casadi::dyn_in (casadi_int ind) |
Get simulator input scheme name by index. More... | |
CASADI_EXPORT std::string | casadi::dyn_out (casadi_int ind) |
Get output scheme name by index. More... | |
CASADI_EXPORT casadi_int | casadi::dyn_n_in () |
Get the number of simulator inputs. More... | |
CASADI_EXPORT casadi_int | casadi::dyn_n_out () |
Get the number of simulator outputs. More... | |
Create an ODE/DAE integrator Solves an initial value problem (IVP) coupled to a terminal value problem with differential equation given as an implicit ODE coupled to an algebraic equation and a set of quadratures:
Initial conditions at t=t0 x(t0) = x0 q(t0) = 0 Forward integration from t=t0 to t=tf der(x) = function(x, z, p, t) Forward ODE 0 = fz(x, z, p, t) Forward algebraic equations der(q) = fq(x, z, p, t) Forward quadratures Terminal conditions at t=tf rx(tf) = rx0 rq(tf) = 0 Backward integration from t=tf to t=t0 der(rx) = gx(rx, rz, rp, x, z, p, t) Backward ODE 0 = gz(rx, rz, rp, x, z, p, t) Backward algebraic equations der(rq) = gq(rx, rz, rp, x, z, p, t) Backward quadratures where we assume that both the forward and backwards integrations are index-1 (i.e. dfz/dz, dgz/drz are invertible) and furthermore that gx, gz and gq have a linear dependency on rx, rz and rp.
Id | Type | Description | Used in |
---|---|---|---|
ad_weight | OT_DOUBLE | Weighting factor for derivative calculation.When there is an option of either using forward or reverse mode directional derivatives, the condition ad_weight*nf<=(1-ad_weight)*na is used where nf and na are estimates of the number of forward/reverse mode directional derivatives needed. By default, ad_weight is calculated automatically, but this can be overridden by setting this option. In particular, 0 means forcing forward mode and 1 forcing reverse mode. Leave unset for (class specific) heuristics. | casadi::FunctionInternal |
ad_weight_sp | OT_DOUBLE | Weighting factor for sparsity pattern calculation calculation.Overrides default behavior. Set to 0 and 1 to force forward and reverse mode respectively. Cf. option "ad_weight". When set to -1, sparsity is completely ignored and dense matrices are used. | casadi::FunctionInternal |
always_inline | OT_BOOL | Force inlining. | casadi::FunctionInternal |
augmented_options | OT_DICT | Options to be passed down to the augmented integrator, if one is constructed. | casadi::Integrator |
cache | OT_DICT | Prepopulate the function cache. Default: empty | casadi::FunctionInternal |
common_options | OT_DICT | Options for auto-generated functions | casadi::OracleFunction |
compiler | OT_STRING | Just-in-time compiler plugin to be used. | casadi::FunctionInternal |
custom_jacobian | OT_FUNCTION | Override CasADi's AD. Use together with 'jac_penalty': 0. Note: Highly experimental. Syntax may break often. | casadi::FunctionInternal |
der_options | OT_DICT | Default options to be used to populate forward_options, reverse_options, and jacobian_options before those options are merged in. | casadi::FunctionInternal |
derivative_of | OT_FUNCTION | The function is a derivative of another function. The type of derivative (directional derivative, Jacobian) is inferred from the function name. | casadi::FunctionInternal |
dump | OT_BOOL | Dump function to file upon first evaluation. [false] | casadi::FunctionInternal |
dump_dir | OT_STRING | Directory to dump inputs/outputs to. Make sure the directory exists [.] | casadi::FunctionInternal |
dump_format | OT_STRING | Choose file format to dump matrices. See DM.from_file [mtx] | casadi::FunctionInternal |
dump_in | OT_BOOL | Dump numerical values of inputs to file (readable with DM.from_file) [default: false] | casadi::FunctionInternal |
dump_out | OT_BOOL | Dump numerical values of outputs to file (readable with DM.from_file) [default: false] | casadi::FunctionInternal |
enable_fd | OT_BOOL | Enable derivative calculation by finite differencing. [default: false]] | casadi::FunctionInternal |
enable_forward | OT_BOOL | Enable derivative calculation using generated functions for Jacobian-times-vector products - typically using forward mode AD - if available. [default: true] | casadi::FunctionInternal |
enable_jacobian | OT_BOOL | Enable derivative calculation using generated functions for Jacobians of all differentiable outputs with respect to all differentiable inputs - if available. [default: true] | casadi::FunctionInternal |
enable_reverse | OT_BOOL | Enable derivative calculation using generated functions for transposed Jacobian-times-vector products - typically using reverse mode AD - if available. [default: true] | casadi::FunctionInternal |
error_on_fail | OT_BOOL | Throw exceptions when function evaluation fails (default true). | casadi::ProtoFunction |
expand | OT_BOOL | Replace MX with SX expressions in problem formulation [false] | casadi::Integrator |
external_transform | OT_VECTORVECTOR | List of external_transform instruction arguments. Default: empty | casadi::FunctionInternal |
fd_method | OT_STRING | Method for finite differencing [default 'central'] | casadi::FunctionInternal |
fd_options | OT_DICT | Options to be passed to the finite difference instance | casadi::FunctionInternal |
forward_options | OT_DICT | Options to be passed to a forward mode constructor | casadi::FunctionInternal |
gather_stats | OT_BOOL | Deprecated option (ignored): Statistics are now always collected. | casadi::FunctionInternal |
grid | OT_DOUBLEVECTOR | [DEPRECATED] Time grid | casadi::Integrator |
input_scheme | OT_STRINGVECTOR | Deprecated option (ignored) | casadi::FunctionInternal |
inputs_check | OT_BOOL | Throw exceptions when the numerical values of the inputs don't make sense | casadi::FunctionInternal |
is_diff_in | OT_BOOLVECTOR | Indicate for each input if it should be differentiable. | casadi::FunctionInternal |
is_diff_out | OT_BOOLVECTOR | Indicate for each output if it should be differentiable. | casadi::FunctionInternal |
jac_penalty | OT_DOUBLE | When requested for a number of forward/reverse directions, it may be cheaper to compute first the full jacobian and then multiply with seeds, rather than obtain the requested directions in a straightforward manner. Casadi uses a heuristic to decide which is cheaper. A high value of 'jac_penalty' makes it less likely for the heurstic to chose the full Jacobian strategy. The special value -1 indicates never to use the full Jacobian strategy | casadi::FunctionInternal |
jacobian_options | OT_DICT | Options to be passed to a Jacobian constructor | casadi::FunctionInternal |
jit | OT_BOOL | Use just-in-time compiler to speed up the evaluation | casadi::FunctionInternal |
jit_cleanup | OT_BOOL | Cleanup up the temporary source file that jit creates. Default: true | casadi::FunctionInternal |
jit_name | OT_STRING | The file name used to write out code. The actual file names used depend on 'jit_temp_suffix' and include extensions. Default: 'jit_tmp' | casadi::FunctionInternal |
jit_options | OT_DICT | Options to be passed to the jit compiler. | casadi::FunctionInternal |
jit_serialize | OT_STRING | Specify behaviour when serializing a jitted function: SOURCE|link|embed. | casadi::FunctionInternal |
jit_temp_suffix | OT_BOOL | Use a temporary (seemingly random) filename suffix for generated code and libraries. This is desired for thread-safety. This behaviour may defeat caching compiler wrappers. Default: true | casadi::FunctionInternal |
max_io | OT_INT | Acceptable number of inputs and outputs. Warn if exceeded. | casadi::FunctionInternal |
max_num_dir | OT_INT | Specify the maximum number of directions for derivative functions. Overrules the builtin optimized_num_dir. | casadi::FunctionInternal |
monitor | OT_STRINGVECTOR | Set of user problem functions to be monitored | casadi::OracleFunction |
nadj | OT_INT | Number of adjoint sensitivities to be calculated [0] | casadi::Integrator |
never_inline | OT_BOOL | Forbid inlining. | casadi::FunctionInternal |
nfwd | OT_INT | Number of forward sensitivities to be calculated [0] | casadi::Integrator |
number_of_finite_elements | OT_INT | Target number of finite elements. The actual number may be higher to accommodate all output times | casadi::Integrator |
output_scheme | OT_STRINGVECTOR | Deprecated option (ignored) | casadi::FunctionInternal |
output_t0 | OT_BOOL | [DEPRECATED] Output the state at the initial time | casadi::Integrator |
post_expand | OT_BOOL | After construction, expand this Function. Default: False | casadi::FunctionInternal |
post_expand_options | OT_DICT | Options to be passed to post-construction expansion. Default: empty | casadi::FunctionInternal |
print_in | OT_BOOL | Print numerical values of inputs [default: false] | casadi::FunctionInternal |
print_out | OT_BOOL | Print numerical values of outputs [default: false] | casadi::FunctionInternal |
print_stats | OT_BOOL | Print out statistics after integration | casadi::Integrator |
print_time | OT_BOOL | print information about execution time. Implies record_time. | casadi::ProtoFunction |
record_time | OT_BOOL | record information about execution time, for retrieval with stats(). | casadi::ProtoFunction |
regularity_check | OT_BOOL | Throw exceptions when NaN or Inf appears during evaluation | casadi::ProtoFunction |
reverse_options | OT_DICT | Options to be passed to a reverse mode constructor | casadi::FunctionInternal |
rootfinder | OT_STRING | An implicit function solver | casadi::Integrator |
rootfinder_options | OT_DICT | Options to be passed to the NLP Solver | casadi::Integrator |
show_eval_warnings | OT_BOOL | Show warnings generated from function evaluations [true] | casadi::OracleFunction |
simplify | OT_BOOL | Implement as MX Function (codegeneratable/serializable) default: false | casadi::Integrator |
simplify_options | OT_DICT | Any options to pass to simplified form Function constructor | casadi::Integrator |
specific_options | OT_DICT | Options for specific auto-generated functions, overwriting the defaults from common_options. Nested dictionary. | casadi::OracleFunction |
t0 | OT_DOUBLE | [DEPRECATED] Beginning of the time horizon | casadi::Integrator |
tf | OT_DOUBLE | [DEPRECATED] End of the time horizon | casadi::Integrator |
user_data | OT_VOIDPTR | A user-defined field that can be used to identify the function or pass additional information | casadi::FunctionInternal |
verbose | OT_BOOL | Verbose evaluation – for debugging | casadi::ProtoFunction |
Full name | Short | Description |
---|---|---|
INTEGRATOR_X0 | x0 | Differential state at the initial time. |
INTEGRATOR_Z0 | z0 | Initial guess for the algebraic variable at the initial time. |
INTEGRATOR_P | p | Parameters. |
INTEGRATOR_U | u | Piecewise constant control, a new control interval starts at each output time. |
INTEGRATOR_ADJ_XF | adj_xf | Adjoint seeds corresponding to the states at the output times. |
INTEGRATOR_ADJ_ZF | adj_zf | Adjoint seeds corresponding to the algebraic variables at the output times. |
INTEGRATOR_ADJ_QF | adj_qf | Adjoint seeds corresponding to the quadratures at the output times. |
Full name | Short | Description |
---|---|---|
INTEGRATOR_XF | xf | Differential state at all output times. |
INTEGRATOR_ZF | zf | Algebraic variable at all output times. |
INTEGRATOR_QF | qf | Quadrature state at all output times. |
INTEGRATOR_ADJ_X0 | adj_x0 | Adjoint sensitivities corresponding to the initial state. |
INTEGRATOR_ADJ_Z0 | adj_z0 | Adjoint sensitivities corresponding to the algebraic variable guess. |
INTEGRATOR_ADJ_P | adj_p | Adjoint sensitivities corresponding to the parameter vector. |
INTEGRATOR_ADJ_U | adj_u | Adjoint sensitivities corresponding to the control vector. |
- cvodes
- idas
- rk
Note: some of the plugins in this list might not be available on your system. Also, there might be extra plugins available to you that are not listed here. You can obtain their documentation with Integrator.doc("myextraplugin")
Interface to CVodes from the Sundials suite.
A call to evaluate will integrate to the end.
You can retrieve the entire state trajectory as follows, after the evaluate call: Call reset. Then call integrate(t_i) and getOuput for a series of times t_i.
Extra doc: https://github.com/casadi/casadi/wiki/L_228
Id | Type | Description |
---|---|---|
abstol | OT_DOUBLE | Absolute tolerence for the IVP solution |
always_recalculate_jacobian | OT_BOOL | Recalculate Jacobian before factorizations, even if Jacobian is current [default: true] |
disable_internal_warnings | OT_BOOL | Disable SUNDIALS internal warning messages |
fsens_all_at_once | OT_BOOL | Calculate all right hand sides of the sensitivity equations at once |
fsens_err_con | OT_BOOL | include the forward sensitivities in all error controls |
interpolation_type | OT_STRING | Type of interpolation for the adjoint sensitivities |
linear_multistep_method | OT_STRING | Integrator scheme: BDF|adams |
linear_solver | OT_STRING | A custom linear solver creator function [default: qr] |
linear_solver_options | OT_DICT | Options to be passed to the linear solver |
max_krylov | OT_INT | Maximum Krylov subspace size |
max_multistep_order | OT_INT | Maximum order for the (variable-order) multistep method |
max_num_steps | OT_INT | Maximum number of integrator steps |
max_order | OT_DOUBLE | Maximum order |
max_step_size | OT_DOUBLE | Max step size [default: 0/inf] |
min_step_size | OT_DOUBLE | Min step size [default: 0/0.0] |
newton_scheme | OT_STRING | Linear solver scheme in the Newton method: DIRECT|gmres|bcgstab|tfqmr |
nonlin_conv_coeff | OT_DOUBLE | Coefficient in the nonlinear convergence test |
nonlinear_solver_iteration | OT_STRING | Nonlinear solver type: NEWTON|functional |
quad_err_con | OT_BOOL | Should the quadratures affect the step size control |
reltol | OT_DOUBLE | Relative tolerence for the IVP solution |
scale_abstol | OT_BOOL | Scale absolute tolerance by nominal value |
second_order_correction | OT_BOOL | Second order correction in the augmented system Jacobian [true] |
sensitivity_method | OT_STRING | Sensitivity method: SIMULTANEOUS|staggered |
step0 | OT_DOUBLE | initial step size [default: 0/estimated] |
steps_per_checkpoint | OT_INT | Number of steps between two consecutive checkpoints |
stop_at_end | OT_BOOL | [DEPRECATED] Stop the integrator at the end of the interval |
use_preconditioner | OT_BOOL | Precondition the iterative solver [default: true] |
Interface to IDAS from the Sundials suite.
Extra doc: https://github.com/casadi/casadi/wiki/L_225
Id | Type | Description |
---|---|---|
abstol | OT_DOUBLE | Absolute tolerence for the IVP solution |
abstolv | OT_DOUBLEVECTOR | Absolute tolerarance for each component |
calc_ic | OT_BOOL | Use IDACalcIC to get consistent initial conditions. |
calc_icB | OT_BOOL | Use IDACalcIC to get consistent initial conditions for backwards system [default: equal to calc_ic]. |
cj_scaling | OT_BOOL | IDAS scaling on cj for the user-defined linear solver module |
constraints | OT_INTVECTOR | Constrain the solution y=[x,z]. 0 (default): no constraint on yi, 1: yi >= 0.0, -1: yi <= 0.0, 2: yi > 0.0, -2: yi < 0.0. |
disable_internal_warnings | OT_BOOL | Disable SUNDIALS internal warning messages |
first_time | OT_DOUBLE | First requested time as a fraction of the time interval |
fsens_err_con | OT_BOOL | include the forward sensitivities in all error controls |
init_xdot | OT_DOUBLEVECTOR | Initial values for the state derivatives |
interpolation_type | OT_STRING | Type of interpolation for the adjoint sensitivities |
linear_solver | OT_STRING | A custom linear solver creator function [default: qr] |
linear_solver_options | OT_DICT | Options to be passed to the linear solver |
max_krylov | OT_INT | Maximum Krylov subspace size |
max_multistep_order | OT_INT | Maximum order for the (variable-order) multistep method |
max_num_steps | OT_INT | Maximum number of integrator steps |
max_order | OT_DOUBLE | Maximum order |
max_step_size | OT_DOUBLE | Maximim step size |
newton_scheme | OT_STRING | Linear solver scheme in the Newton method: DIRECT|gmres|bcgstab|tfqmr |
nonlin_conv_coeff | OT_DOUBLE | Coefficient in the nonlinear convergence test |
quad_err_con | OT_BOOL | Should the quadratures affect the step size control |
reltol | OT_DOUBLE | Relative tolerence for the IVP solution |
scale_abstol | OT_BOOL | Scale absolute tolerance by nominal value |
second_order_correction | OT_BOOL | Second order correction in the augmented system Jacobian [true] |
sensitivity_method | OT_STRING | Sensitivity method: SIMULTANEOUS|staggered |
step0 | OT_DOUBLE | initial step size [default: 0/estimated] |
steps_per_checkpoint | OT_INT | Number of steps between two consecutive checkpoints |
stop_at_end | OT_BOOL | [DEPRECATED] Stop the integrator at the end of the interval |
suppress_algebraic | OT_BOOL | Suppress algebraic variables in the error testing |
use_preconditioner | OT_BOOL | Precondition the iterative solver [default: true] |
Fixed-step implicit Runge-Kutta integrator ODE/DAE integrator based on collocation schemes
The method is still under development
Extra doc: https://github.com/casadi/casadi/wiki/L_234
Id | Type | Description |
---|---|---|
collocation_scheme | OT_STRING | Collocation scheme: radau|legendre |
interpolation_order | OT_INT | Order of the interpolating polynomials |
number_of_finite_elements | OT_INT | Target number of finite elements. The actual number may be higher to accommodate all output times |
rootfinder | OT_STRING | An implicit function solver |
rootfinder_options | OT_DICT | Options to be passed to the NLP Solver |
simplify | OT_BOOL | Implement as MX Function (codegeneratable/serializable) default: false |
simplify_options | OT_DICT | Any options to pass to simplified form Function constructor |
Fixed-step explicit Runge-Kutta integrator for ODEs Currently implements RK4.
The method is still under development
Extra doc: https://github.com/casadi/casadi/wiki/L_23a
Extra doc: https://github.com/casadi/casadi/wiki/L_21k
CASADI_EXPORT std::string casadi::doc_integrator | ( | const std::string & | name | ) |
CASADI_EXPORT std::vector<std::string> casadi::dyn_in | ( | ) |
Extra doc: https://github.com/casadi/casadi/wiki/L_25p
CASADI_EXPORT std::string casadi::dyn_in | ( | casadi_int | ind | ) |
Extra doc: https://github.com/casadi/casadi/wiki/L_25r
CASADI_EXPORT casadi_int casadi::dyn_n_in | ( | ) |
Extra doc: https://github.com/casadi/casadi/wiki/L_25t
CASADI_EXPORT casadi_int casadi::dyn_n_out | ( | ) |
Extra doc: https://github.com/casadi/casadi/wiki/L_25u
CASADI_EXPORT std::vector<std::string> casadi::dyn_out | ( | ) |
Extra doc: https://github.com/casadi/casadi/wiki/L_25q
CASADI_EXPORT std::string casadi::dyn_out | ( | casadi_int | ind | ) |
Extra doc: https://github.com/casadi/casadi/wiki/L_25s
CASADI_EXPORT bool casadi::has_integrator | ( | const std::string & | name | ) |
CASADI_EXPORT Function casadi::integrator | ( | const std::string & | name, |
const std::string & | solver, | ||
const Function & | dae, | ||
const Dict & | opts = Dict() |
||
) |
CASADI_EXPORT Function casadi::integrator | ( | const std::string & | name, |
const std::string & | solver, | ||
const Function & | dae, | ||
double | t0, | ||
const std::vector< double > & | tout, | ||
const Dict & | opts = Dict() |
||
) |
CASADI_EXPORT Function casadi::integrator | ( | const std::string & | name, |
const std::string & | solver, | ||
const Function & | dae, | ||
double | t0, | ||
double | tf, | ||
const Dict & | opts = Dict() |
||
) |
CASADI_EXPORT Function casadi::integrator | ( | const std::string & | name, |
const std::string & | solver, | ||
const MXDict & | dae, | ||
const Dict & | opts = Dict() |
||
) |
CASADI_EXPORT Function casadi::integrator | ( | const std::string & | name, |
const std::string & | solver, | ||
const MXDict & | dae, | ||
double | t0, | ||
const std::vector< double > & | tout, | ||
const Dict & | opts = Dict() |
||
) |
CASADI_EXPORT Function casadi::integrator | ( | const std::string & | name, |
const std::string & | solver, | ||
const MXDict & | dae, | ||
double | t0, | ||
double | tf, | ||
const Dict & | opts = Dict() |
||
) |
CASADI_EXPORT Function casadi::integrator | ( | const std::string & | name, |
const std::string & | solver, | ||
const SXDict & | dae, | ||
const Dict & | opts = Dict() |
||
) |
CASADI_EXPORT Function casadi::integrator | ( | const std::string & | name, |
const std::string & | solver, | ||
const SXDict & | dae, | ||
double | t0, | ||
const std::vector< double > & | tout, | ||
const Dict & | opts = Dict() |
||
) |
CASADI_EXPORT Function casadi::integrator | ( | const std::string & | name, |
const std::string & | solver, | ||
const SXDict & | dae, | ||
double | t0, | ||
double | tf, | ||
const Dict & | opts = Dict() |
||
) |
CASADI_EXPORT std::vector<std::string> casadi::integrator_in | ( | ) |
Extra doc: https://github.com/casadi/casadi/wiki/L_7b
CASADI_EXPORT std::string casadi::integrator_in | ( | casadi_int | ind | ) |
Extra doc: https://github.com/casadi/casadi/wiki/L_7d
CASADI_EXPORT casadi_int casadi::integrator_n_in | ( | ) |
Extra doc: https://github.com/casadi/casadi/wiki/L_7f
CASADI_EXPORT casadi_int casadi::integrator_n_out | ( | ) |
Extra doc: https://github.com/casadi/casadi/wiki/L_7g
CASADI_EXPORT std::vector<std::string> casadi::integrator_out | ( | ) |
Extra doc: https://github.com/casadi/casadi/wiki/L_7c
CASADI_EXPORT std::string casadi::integrator_out | ( | casadi_int | ind | ) |
Extra doc: https://github.com/casadi/casadi/wiki/L_7e
CASADI_EXPORT void casadi::load_integrator | ( | const std::string & | name | ) |