Functions
Title

Functions

Function casadi::integrator (const std::string &name, const std::string &solver, const SXDict &dae, const Dict &opts)
 
Function casadi::integrator (const std::string &name, const std::string &solver, const MXDict &dae, const Dict &opts)
 
Function casadi::integrator (const std::string &name, const std::string &solver, const Function &dae, const Dict &opts)
 
Function casadi::integrator (const std::string &name, const std::string &solver, const SXDict &dae, double t0, const std::vector< double > &tout, const Dict &opts)
 
Function casadi::integrator (const std::string &name, const std::string &solver, const MXDict &dae, double t0, const std::vector< double > &tout, const Dict &opts)
 
Function casadi::integrator (const std::string &name, const std::string &solver, const Function &dae, double t0, const std::vector< double > &tout, const Dict &opts)
 
Function casadi::integrator (const std::string &name, const std::string &solver, const SXDict &dae, double t0, double tf, const Dict &opts)
 
Function casadi::integrator (const std::string &name, const std::string &solver, const MXDict &dae, double t0, double tf, const Dict &opts)
 
Function casadi::integrator (const std::string &name, const std::string &solver, const Function &dae, double t0, double tf, const Dict &opts)
 
bool casadi::has_integrator (const std::string &name)
 Check if a particular plugin is available. More...
 
void casadi::load_integrator (const std::string &name)
 Explicitly load a plugin dynamically. More...
 
std::string casadi::doc_integrator (const std::string &name)
 Get the documentation string for a plugin. More...
 
std::vector< std::string > casadi::integrator_in ()
 Get input scheme of integrators. More...
 
std::vector< std::string > casadi::integrator_out ()
 Get integrator output scheme of integrators. More...
 
std::string casadi::integrator_in (casadi_int ind)
 Get integrator input scheme name by index. More...
 
std::string casadi::integrator_out (casadi_int ind)
 Get output scheme name by index. More...
 
casadi_int casadi::integrator_n_in ()
 Get the number of integrator inputs. More...
 
casadi_int casadi::integrator_n_out ()
 Get the number of integrator outputs. More...
 
std::vector< std::string > casadi::dyn_in ()
 Get input scheme of a DAE function. More...
 
std::vector< std::string > casadi::dyn_out ()
 Get output scheme of a DAE function. More...
 
std::string casadi::dyn_in (casadi_int ind)
 Get input scheme of a DAE function by index. More...
 
std::string casadi::dyn_out (casadi_int ind)
 Get output scheme of a DAE function by index. More...
 
casadi_int casadi::dyn_n_in ()
 Get the number of inputs for a DAE function. More...
 
casadi_int casadi::dyn_n_out ()
 Get the number of outputs for a DAE function. More...
 
std::vector< std::string > casadi::event_in ()
 Get input scheme of an event transition function. More...
 
std::vector< std::string > casadi::event_out ()
 Get output scheme of an event transition functions. More...
 

Detailed Description

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.

General information


List of available options
IdTypeDescriptionUsed in
ad_weightOT_DOUBLEWeighting 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_spOT_DOUBLEWeighting 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_inlineOT_BOOLForce inlining.casadi::FunctionInternal
augmented_optionsOT_DICTOptions to be passed down to the augmented integrator, if one is constructedcasadi::Integrator
cacheOT_DICTPrepopulate the function cache. Default: emptycasadi::FunctionInternal
common_optionsOT_DICTOptions for auto-generated functionscasadi::OracleFunction
compilerOT_STRINGJust-in-time compiler plugin to be used.casadi::FunctionInternal
custom_jacobianOT_FUNCTIONOverride CasADi's AD. Use together with 'jac_penalty': 0. Note: Highly experimental. Syntax may break often.casadi::FunctionInternal
der_optionsOT_DICTDefault options to be used to populate forward_options, reverse_options, and jacobian_options before those options are merged in.casadi::FunctionInternal
derivative_ofOT_FUNCTIONThe function is a derivative of another function. The type of derivative (directional derivative, Jacobian) is inferred from the function name.casadi::FunctionInternal
dumpOT_BOOLDump function to file upon first evaluation. [false]casadi::FunctionInternal
dump_dirOT_STRINGDirectory to dump inputs/outputs to. Make sure the directory exists [.]casadi::FunctionInternal
dump_formatOT_STRINGChoose file format to dump matrices. See DM.from_file [mtx]casadi::FunctionInternal
dump_inOT_BOOLDump numerical values of inputs to file (readable with DM.from_file) [default: false]casadi::FunctionInternal
dump_outOT_BOOLDump numerical values of outputs to file (readable with DM.from_file) [default: false]casadi::FunctionInternal
enable_fdOT_BOOLEnable derivative calculation by finite differencing. [default: false]]casadi::FunctionInternal
enable_forwardOT_BOOLEnable derivative calculation using generated functions for Jacobian-times-vector products - typically using forward mode AD - if available. [default: true]casadi::FunctionInternal
enable_jacobianOT_BOOLEnable derivative calculation using generated functions for Jacobians of all differentiable outputs with respect to all differentiable inputs - if available. [default: true]casadi::FunctionInternal
enable_reverseOT_BOOLEnable derivative calculation using generated functions for transposed Jacobian-times-vector products - typically using reverse mode AD - if available. [default: true]casadi::FunctionInternal
error_on_failOT_BOOLThrow exceptions when function evaluation fails (default true).casadi::ProtoFunction
event_tolOT_DOUBLETermination tolerance for the event iterationcasadi::Integrator
expandOT_BOOLReplace MX with SX expressions in problem formulation [false] This happens before creating derivatives unless indicated by postpone_expandcasadi::OracleFunction
external_transformOT_VECTORVECTORList of external_transform instruction arguments. Default: emptycasadi::FunctionInternal
fd_methodOT_STRINGMethod for finite differencing [default 'central']casadi::FunctionInternal
fd_optionsOT_DICTOptions to be passed to the finite difference instancecasadi::FunctionInternal
forward_optionsOT_DICTOptions to be passed to a forward mode constructorcasadi::FunctionInternal
gather_statsOT_BOOLDeprecated option (ignored): Statistics are now always collected.casadi::FunctionInternal
gridOT_DOUBLEVECTOR[DEPRECATED] Time gridcasadi::Integrator
input_schemeOT_STRINGVECTORDeprecated option (ignored)casadi::FunctionInternal
inputs_checkOT_BOOLThrow exceptions when the numerical values of the inputs don't make sensecasadi::FunctionInternal
is_diff_inOT_BOOLVECTORIndicate for each input if it should be differentiable.casadi::FunctionInternal
is_diff_outOT_BOOLVECTORIndicate for each output if it should be differentiable.casadi::FunctionInternal
jac_penaltyOT_DOUBLEWhen 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 strategycasadi::FunctionInternal
jacobian_optionsOT_DICTOptions to be passed to a Jacobian constructorcasadi::FunctionInternal
jitOT_BOOLUse just-in-time compiler to speed up the evaluationcasadi::FunctionInternal
jit_cleanupOT_BOOLCleanup up the temporary source file that jit creates. Default: truecasadi::FunctionInternal
jit_nameOT_STRINGThe 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_optionsOT_DICTOptions to be passed to the jit compiler.casadi::FunctionInternal
jit_serializeOT_STRINGSpecify behaviour when serializing a jitted function: SOURCE|link|embed.casadi::FunctionInternal
jit_temp_suffixOT_BOOLUse 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: truecasadi::FunctionInternal
max_event_iterOT_INTMaximum number of iterations to zero in on a single eventcasadi::Integrator
max_eventsOT_INTMaximum total number of eventscasadi::Integrator
max_ioOT_INTAcceptable number of inputs and outputs. Warn if exceeded.casadi::FunctionInternal
max_num_dirOT_INTSpecify the maximum number of directions for derivative functions. Overrules the builtin optimized_num_dir.casadi::FunctionInternal
monitorOT_STRINGVECTORSet of user problem functions to be monitoredcasadi::OracleFunction
nadjOT_INTNumber of adjoint sensitivities to be calculated [0]casadi::Integrator
never_inlineOT_BOOLForbid inlining.casadi::FunctionInternal
nfwdOT_INTNumber of forward sensitivities to be calculated [0]casadi::Integrator
number_of_finite_elementsOT_INTTarget number of finite elements. The actual number may be higher to accommodate all output timescasadi::Integrator
output_schemeOT_STRINGVECTORDeprecated option (ignored)casadi::FunctionInternal
output_t0OT_BOOL[DEPRECATED] Output the state at the initial timecasadi::Integrator
post_expandOT_BOOLAfter construction, expand this Function. Default: Falsecasadi::FunctionInternal
post_expand_optionsOT_DICTOptions to be passed to post-construction expansion. Default: emptycasadi::FunctionInternal
postpone_expandOT_BOOLWhen expand is active, postpone it until after creation of derivatives. Default: Falsecasadi::OracleFunction
print_inOT_BOOLPrint numerical values of inputs [default: false]casadi::FunctionInternal
print_outOT_BOOLPrint numerical values of outputs [default: false]casadi::FunctionInternal
print_statsOT_BOOLPrint out statistics after integrationcasadi::Integrator
print_timeOT_BOOLprint information about execution time. Implies record_time.casadi::ProtoFunction
record_timeOT_BOOLrecord information about execution time, for retrieval with stats().casadi::ProtoFunction
regularity_checkOT_BOOLThrow exceptions when NaN or Inf appears during evaluationcasadi::ProtoFunction
reverse_optionsOT_DICTOptions to be passed to a reverse mode constructorcasadi::FunctionInternal
rootfinderOT_STRINGAn implicit function solvercasadi::Integrator
rootfinder_optionsOT_DICTOptions to be passed to the NLP Solvercasadi::Integrator
show_eval_warningsOT_BOOLShow warnings generated from function evaluations [true]casadi::OracleFunction
simplifyOT_BOOLImplement as MX Function (codegeneratable/serializable) default: falsecasadi::Integrator
simplify_optionsOT_DICTAny options to pass to simplified form Function constructorcasadi::Integrator
specific_optionsOT_DICTOptions for specific auto-generated functions, overwriting the defaults from common_options. Nested dictionary.casadi::OracleFunction
t0OT_DOUBLE[DEPRECATED] Beginning of the time horizoncasadi::Integrator
tfOT_DOUBLE[DEPRECATED] End of the time horizoncasadi::Integrator
transitionOT_FUNCTIONFunction to be called a zero-crossing eventscasadi::Integrator
user_dataOT_VOIDPTRA user-defined field that can be used to identify the function or pass additional informationcasadi::FunctionInternal
verboseOT_BOOLVerbose evaluation – for debuggingcasadi::ProtoFunction

Input scheme: casadi::IntegratorInput (INTEGRATOR_NUM_IN = 7)
Full nameShortDescription
INTEGRATOR_X0x0Differential state at the initial time.
INTEGRATOR_Z0z0Initial guess for the algebraic variable at the initial time.
INTEGRATOR_PpParameters.
INTEGRATOR_UuPiecewise constant control, a new control interval starts at each output time.
INTEGRATOR_ADJ_XFadj_xfAdjoint seeds corresponding to the states at the output times.
INTEGRATOR_ADJ_ZFadj_zfAdjoint seeds corresponding to the algebraic variables at the output times.
INTEGRATOR_ADJ_QFadj_qfAdjoint seeds corresponding to the quadratures at the output times.

Output scheme: casadi::IntegratorOutput (INTEGRATOR_NUM_OUT = 7)
Full nameShortDescription
INTEGRATOR_XFxfDifferential state at all output times.
INTEGRATOR_ZFzfAlgebraic variable at all output times.
INTEGRATOR_QFqfQuadrature state at all output times.
INTEGRATOR_ADJ_X0adj_x0Adjoint sensitivities corresponding to the initial state.
INTEGRATOR_ADJ_Z0adj_z0Adjoint sensitivities corresponding to the algebraic variable guess.
INTEGRATOR_ADJ_Padj_pAdjoint sensitivities corresponding to the parameter vector.
INTEGRATOR_ADJ_Uadj_uAdjoint sensitivities corresponding to the control vector.

List of plugins

- cvodes

- idas

- collocation

- 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")


cvodes

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


List of available options
IdTypeDescription
abstolOT_DOUBLEAbsolute tolerence for the IVP solution
always_recalculate_jacobianOT_BOOLRecalculate Jacobian before factorizations, even if Jacobian is current [default: true]
disable_internal_warningsOT_BOOLDisable SUNDIALS internal warning messages
fsens_all_at_onceOT_BOOLCalculate all right hand sides of the sensitivity equations at once
fsens_err_conOT_BOOLinclude the forward sensitivities in all error controls
interpolation_typeOT_STRINGType of interpolation for the adjoint sensitivities
linear_multistep_methodOT_STRINGIntegrator scheme: BDF|adams
linear_solverOT_STRINGA custom linear solver creator function [default: qr]
linear_solver_optionsOT_DICTOptions to be passed to the linear solver
max_krylovOT_INTMaximum Krylov subspace size
max_multistep_orderOT_INTMaximum order for the (variable-order) multistep method
max_num_stepsOT_INTMaximum number of integrator steps
max_orderOT_DOUBLEMaximum order
max_step_sizeOT_DOUBLEMax step size [default: 0/inf]
min_step_sizeOT_DOUBLEMin step size [default: 0/0.0]
newton_schemeOT_STRINGLinear solver scheme in the Newton method: DIRECT|gmres|bcgstab|tfqmr
nonlin_conv_coeffOT_DOUBLECoefficient in the nonlinear convergence test
nonlinear_solver_iterationOT_STRINGNonlinear solver type: NEWTON|functional
quad_err_conOT_BOOLShould the quadratures affect the step size control
reltolOT_DOUBLERelative tolerence for the IVP solution
scale_abstolOT_BOOLScale absolute tolerance by nominal value
second_order_correctionOT_BOOLSecond order correction in the augmented system Jacobian [true]
sensitivity_methodOT_STRINGSensitivity method: SIMULTANEOUS|staggered
step0OT_DOUBLEinitial step size [default: 0/estimated]
steps_per_checkpointOT_INTNumber of steps between two consecutive checkpoints
stop_at_endOT_BOOL[DEPRECATED] Stop the integrator at the end of the interval
use_preconditionerOT_BOOLPrecondition the iterative solver [default: true]

idas

Interface to IDAS from the Sundials suite.

Extra doc: https://github.com/casadi/casadi/wiki/L_225


List of available options
IdTypeDescription
abstolOT_DOUBLEAbsolute tolerence for the IVP solution
abstolvOT_DOUBLEVECTORAbsolute tolerarance for each component
calc_icOT_BOOLUse IDACalcIC to get consistent initial conditions.
calc_icBOT_BOOLUse IDACalcIC to get consistent initial conditions for backwards system [default: equal to calc_ic].
cj_scalingOT_BOOLIDAS scaling on cj for the user-defined linear solver module
constraintsOT_INTVECTORConstrain 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_warningsOT_BOOLDisable SUNDIALS internal warning messages
first_timeOT_DOUBLEFirst requested time as a fraction of the time interval
fsens_err_conOT_BOOLinclude the forward sensitivities in all error controls
init_xdotOT_DOUBLEVECTORInitial values for the state derivatives
interpolation_typeOT_STRINGType of interpolation for the adjoint sensitivities
linear_solverOT_STRINGA custom linear solver creator function [default: qr]
linear_solver_optionsOT_DICTOptions to be passed to the linear solver
max_krylovOT_INTMaximum Krylov subspace size
max_multistep_orderOT_INTMaximum order for the (variable-order) multistep method
max_num_stepsOT_INTMaximum number of integrator steps
max_orderOT_DOUBLEMaximum order
max_step_sizeOT_DOUBLEMaximim step size
newton_schemeOT_STRINGLinear solver scheme in the Newton method: DIRECT|gmres|bcgstab|tfqmr
nonlin_conv_coeffOT_DOUBLECoefficient in the nonlinear convergence test
quad_err_conOT_BOOLShould the quadratures affect the step size control
reltolOT_DOUBLERelative tolerence for the IVP solution
scale_abstolOT_BOOLScale absolute tolerance by nominal value
second_order_correctionOT_BOOLSecond order correction in the augmented system Jacobian [true]
sensitivity_methodOT_STRINGSensitivity method: SIMULTANEOUS|staggered
step0OT_DOUBLEinitial step size [default: 0/estimated]
steps_per_checkpointOT_INTNumber of steps between two consecutive checkpoints
stop_at_endOT_BOOL[DEPRECATED] Stop the integrator at the end of the interval
suppress_algebraicOT_BOOLSuppress algebraic variables in the error testing
use_preconditionerOT_BOOLPrecondition the iterative solver [default: true]

collocation

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


List of available options
IdTypeDescription
collocation_schemeOT_STRINGCollocation scheme: radau|legendre
interpolation_orderOT_INTOrder of the interpolating polynomials
number_of_finite_elementsOT_INTTarget number of finite elements. The actual number may be higher to accommodate all output times
rootfinderOT_STRINGAn implicit function solver
rootfinder_optionsOT_DICTOptions to be passed to the NLP Solver
simplifyOT_BOOLImplement as MX Function (codegeneratable/serializable) default: false
simplify_optionsOT_DICTAny options to pass to simplified form Function constructor

rk

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

Author
Joel Andersson
Date
2011-2015

Extra doc: https://github.com/casadi/casadi/wiki/L_21k

Function Documentation

◆ doc_integrator()

CASADI_EXPORT std::string casadi::doc_integrator ( const std::string &  name)

Definition at line 130 of file integrator.cpp.

130  {
131  return Integrator::getPlugin(name).doc;
132 }

References casadi::PluginInterface< Integrator >::getPlugin().

◆ dyn_in() [1/2]

CASADI_EXPORT std::vector< std::string > casadi::dyn_in ( )

◆ dyn_in() [2/2]

CASADI_EXPORT std::string casadi::dyn_in ( casadi_int  ind)

Extra doc: https://github.com/casadi/casadi/wiki/L_25r

Definition at line 240 of file integrator.cpp.

240  {
241  return to_string(static_cast<DynIn>(ind));
242 }
DynIn
Inputs of the symbolic representation of the DAE.
Definition: integrator.hpp:190
std::string to_string(TypeFmi2 v)

References casadi::to_string().

◆ dyn_n_in()

CASADI_EXPORT casadi_int casadi::dyn_n_in ( )

Extra doc: https://github.com/casadi/casadi/wiki/L_25t

Definition at line 248 of file integrator.cpp.

248  {
249  return DYN_NUM_IN;
250 }
@ DYN_NUM_IN
Definition: integrator.hpp:196

References casadi::DYN_NUM_IN.

◆ dyn_n_out()

CASADI_EXPORT casadi_int casadi::dyn_n_out ( )

Extra doc: https://github.com/casadi/casadi/wiki/L_25u

Definition at line 252 of file integrator.cpp.

252  {
253  return DYN_NUM_OUT;
254 }
@ DYN_NUM_OUT
Definition: integrator.hpp:204

References casadi::DYN_NUM_OUT.

◆ dyn_out() [1/2]

CASADI_EXPORT std::vector< std::string > casadi::dyn_out ( )

◆ dyn_out() [2/2]

CASADI_EXPORT std::string casadi::dyn_out ( casadi_int  ind)

Extra doc: https://github.com/casadi/casadi/wiki/L_25s

Definition at line 244 of file integrator.cpp.

244  {
245  return to_string(static_cast<DynOut>(ind));
246 }
DynOut
Outputs of the symbolic representation of the DAE.
Definition: integrator.hpp:199

References casadi::to_string().

◆ event_in()

CASADI_EXPORT std::vector< std::string > casadi::event_in ( )

Extra doc: https://github.com/casadi/casadi/wiki/L_2b4

Definition at line 256 of file integrator.cpp.

256  {
257  return enum_names<EventIn>();
258 }

Referenced by casadi::DaeBuilderInternal::transition().

◆ event_out()

CASADI_EXPORT std::vector< std::string > casadi::event_out ( )

Extra doc: https://github.com/casadi/casadi/wiki/L_2b5

Definition at line 260 of file integrator.cpp.

260  {
261  return enum_names<EventOut>();
262 }

Referenced by casadi::DaeBuilderInternal::transition().

◆ has_integrator()

CASADI_EXPORT bool casadi::has_integrator ( const std::string &  name)

Definition at line 122 of file integrator.cpp.

122  {
123  return Integrator::has_plugin(name);
124 }

References casadi::PluginInterface< Integrator >::has_plugin().

◆ integrator() [1/9]

CASADI_EXPORT Function casadi::integrator ( const std::string &  name,
const std::string &  solver,
const Function dae,
const Dict opts 
)

Definition at line 144 of file integrator.cpp.

145  {
146  return integrator(name, solver, dae, 0.0, std::vector<double>{1.0}, opts);
147 }
Function integrator(const std::string &name, const std::string &solver, const Function &dae, double t0, double tf, const Dict &opts)
Definition: integrator.cpp:179

References casadi::integrator().

◆ integrator() [2/9]

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 
)

Definition at line 159 of file integrator.cpp.

160  {
161  // Make sure that dae is sound
162  if (dae.has_free()) {
163  casadi_error("Cannot create '" + name + "' since " + str(dae.get_free()) + " are free.");
164  }
165  Integrator* intg = Integrator::getPlugin(solver).creator(name, dae, t0, tout);
166  return intg->create_advanced(opts);
167 }
std::string str(const T &v)
String representation, any type.

References casadi::Integrator::create_advanced(), casadi::Function::get_free(), casadi::PluginInterface< Integrator >::getPlugin(), casadi::Function::has_free(), and casadi::str().

◆ integrator() [3/9]

CASADI_EXPORT Function casadi::integrator ( const std::string &  name,
const std::string &  solver,
const Function dae,
double  t0,
double  tf,
const Dict opts 
)

Definition at line 179 of file integrator.cpp.

180  {
181  return integrator(name, solver, dae, t0, std::vector<double>{tf}, opts);
182 }

References casadi::integrator().

◆ integrator() [4/9]

CASADI_EXPORT Function casadi::integrator ( const std::string &  name,
const std::string &  solver,
const MXDict dae,
const Dict opts 
)

Definition at line 139 of file integrator.cpp.

140  {
141  return integrator(name, solver, dae, 0.0, std::vector<double>{1.0}, opts);
142 }

References casadi::integrator().

◆ integrator() [5/9]

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 
)

Definition at line 154 of file integrator.cpp.

155  {
156  return integrator(name, solver, Integrator::map2oracle("dae", dae), t0, tout, opts);
157 }

References casadi::integrator(), and casadi::Integrator::map2oracle().

◆ integrator() [6/9]

CASADI_EXPORT Function casadi::integrator ( const std::string &  name,
const std::string &  solver,
const MXDict dae,
double  t0,
double  tf,
const Dict opts 
)

Definition at line 174 of file integrator.cpp.

175  {
176  return integrator(name, solver, dae, t0, std::vector<double>{tf}, opts);
177 }

References casadi::integrator().

◆ integrator() [7/9]

CASADI_EXPORT Function casadi::integrator ( const std::string &  name,
const std::string &  solver,
const SXDict dae,
const Dict opts = Dict() 
)
Examples
integrators/idas.py, and integrators/tolerance.py.

Definition at line 134 of file integrator.cpp.

135  {
136  return integrator(name, solver, dae, 0.0, std::vector<double>{1.0}, opts);
137 }

Referenced by casadi::Integrator::get_forward(), casadi::Integrator::get_reverse(), casadi::integrator(), and casadi::simpleIntegrator().

◆ integrator() [8/9]

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 
)

Definition at line 149 of file integrator.cpp.

150  {
151  return integrator(name, solver, Integrator::map2oracle("dae", dae), t0, tout, opts);
152 }

References casadi::integrator(), and casadi::Integrator::map2oracle().

◆ integrator() [9/9]

CASADI_EXPORT Function casadi::integrator ( const std::string &  name,
const std::string &  solver,
const SXDict dae,
double  t0,
double  tf,
const Dict opts 
)

Definition at line 169 of file integrator.cpp.

170  {
171  return integrator(name, solver, dae, t0, std::vector<double>{tf}, opts);
172 }

References casadi::integrator().

◆ integrator_in() [1/2]

CASADI_EXPORT std::vector< std::string > casadi::integrator_in ( )

Extra doc: https://github.com/casadi/casadi/wiki/L_7b

Definition at line 184 of file integrator.cpp.

184  {
185  std::vector<std::string> ret(integrator_n_in());
186  for (size_t i=0; i<ret.size(); ++i) ret[i]=integrator_in(i);
187  return ret;
188 }
std::string integrator_in(casadi_int ind)
Get integrator input scheme name by index.
Definition: integrator.cpp:196
casadi_int integrator_n_in()
Get the number of integrator inputs.
Definition: integrator.cpp:224

References casadi::integrator_n_in().

Referenced by casadi::FixedStepIntegrator::create_advanced(), casadi::Integrator::get_forward(), casadi::Integrator::get_name_in(), and casadi::Integrator::get_reverse().

◆ integrator_in() [2/2]

CASADI_EXPORT std::string casadi::integrator_in ( casadi_int  ind)

Extra doc: https://github.com/casadi/casadi/wiki/L_7d

Definition at line 196 of file integrator.cpp.

196  {
197  switch (static_cast<IntegratorInput>(ind)) {
198  case INTEGRATOR_X0: return "x0";
199  case INTEGRATOR_Z0: return "z0";
200  case INTEGRATOR_P: return "p";
201  case INTEGRATOR_U: return "u";
202  case INTEGRATOR_ADJ_XF: return "adj_xf";
203  case INTEGRATOR_ADJ_ZF: return "adj_zf";
204  case INTEGRATOR_ADJ_QF: return "adj_qf";
205  case INTEGRATOR_NUM_IN: break;
206  }
207  return std::string();
208 }
IntegratorInput
Input arguments of an integrator.
Definition: integrator.hpp:223
@ INTEGRATOR_U
Piecewise constant control, a new control interval starts at each output time.
Definition: integrator.hpp:231
@ INTEGRATOR_ADJ_QF
Adjoint seeds corresponding to the quadratures at the output times.
Definition: integrator.hpp:237
@ INTEGRATOR_ADJ_ZF
Adjoint seeds corresponding to the algebraic variables at the output times.
Definition: integrator.hpp:235
@ INTEGRATOR_P
Parameters.
Definition: integrator.hpp:229
@ INTEGRATOR_ADJ_XF
Adjoint seeds corresponding to the states at the output times.
Definition: integrator.hpp:233
@ INTEGRATOR_Z0
Initial guess for the algebraic variable at the initial time.
Definition: integrator.hpp:227
@ INTEGRATOR_NUM_IN
Number of input arguments of an integrator.
Definition: integrator.hpp:239
@ INTEGRATOR_X0
Differential state at the initial time.
Definition: integrator.hpp:225

References casadi::INTEGRATOR_ADJ_QF, casadi::INTEGRATOR_ADJ_XF, casadi::INTEGRATOR_ADJ_ZF, casadi::INTEGRATOR_NUM_IN, casadi::INTEGRATOR_P, casadi::INTEGRATOR_U, casadi::INTEGRATOR_X0, and casadi::INTEGRATOR_Z0.

◆ integrator_n_in()

CASADI_EXPORT casadi_int casadi::integrator_n_in ( )

Extra doc: https://github.com/casadi/casadi/wiki/L_7f

Definition at line 224 of file integrator.cpp.

224  {
225  return INTEGRATOR_NUM_IN;
226 }

References casadi::INTEGRATOR_NUM_IN.

Referenced by casadi::integrator_in().

◆ integrator_n_out()

CASADI_EXPORT casadi_int casadi::integrator_n_out ( )

Extra doc: https://github.com/casadi/casadi/wiki/L_7g

Definition at line 228 of file integrator.cpp.

228  {
229  return INTEGRATOR_NUM_OUT;
230 }
@ INTEGRATOR_NUM_OUT
Number of output arguments of an integrator.
Definition: integrator.hpp:259

References casadi::INTEGRATOR_NUM_OUT.

Referenced by casadi::integrator_out().

◆ integrator_out() [1/2]

CASADI_EXPORT std::vector< std::string > casadi::integrator_out ( )

Extra doc: https://github.com/casadi/casadi/wiki/L_7c

Definition at line 190 of file integrator.cpp.

190  {
191  std::vector<std::string> ret(integrator_n_out());
192  for (size_t i=0; i<ret.size(); ++i) ret[i]=integrator_out(i);
193  return ret;
194 }
casadi_int integrator_n_out()
Get the number of integrator outputs.
Definition: integrator.cpp:228
std::string integrator_out(casadi_int ind)
Get output scheme name by index.
Definition: integrator.cpp:210

References casadi::integrator_n_out().

Referenced by casadi::FixedStepIntegrator::create_advanced(), casadi::Integrator::get_forward(), casadi::Integrator::get_name_out(), and casadi::Integrator::get_reverse().

◆ integrator_out() [2/2]

CASADI_EXPORT std::string casadi::integrator_out ( casadi_int  ind)

Extra doc: https://github.com/casadi/casadi/wiki/L_7e

Definition at line 210 of file integrator.cpp.

210  {
211  switch (static_cast<IntegratorOutput>(ind)) {
212  case INTEGRATOR_XF: return "xf";
213  case INTEGRATOR_ZF: return "zf";
214  case INTEGRATOR_QF: return "qf";
215  case INTEGRATOR_ADJ_X0: return "adj_x0";
216  case INTEGRATOR_ADJ_Z0: return "adj_z0";
217  case INTEGRATOR_ADJ_P: return "adj_p";
218  case INTEGRATOR_ADJ_U: return "adj_u";
219  case INTEGRATOR_NUM_OUT: break;
220  }
221  return std::string();
222 }
IntegratorOutput
Output arguments of an integrator.
Definition: integrator.hpp:243
@ INTEGRATOR_ADJ_U
Adjoint sensitivities corresponding to the control vector.
Definition: integrator.hpp:257
@ INTEGRATOR_ADJ_Z0
Adjoint sensitivities corresponding to the algebraic variable guess.
Definition: integrator.hpp:253
@ INTEGRATOR_QF
Quadrature state at all output times.
Definition: integrator.hpp:249
@ INTEGRATOR_ADJ_P
Adjoint sensitivities corresponding to the parameter vector.
Definition: integrator.hpp:255
@ INTEGRATOR_ZF
Algebraic variable at all output times.
Definition: integrator.hpp:247
@ INTEGRATOR_XF
Differential state at all output times.
Definition: integrator.hpp:245
@ INTEGRATOR_ADJ_X0
Adjoint sensitivities corresponding to the initial state.
Definition: integrator.hpp:251

References casadi::INTEGRATOR_ADJ_P, casadi::INTEGRATOR_ADJ_U, casadi::INTEGRATOR_ADJ_X0, casadi::INTEGRATOR_ADJ_Z0, casadi::INTEGRATOR_NUM_OUT, casadi::INTEGRATOR_QF, casadi::INTEGRATOR_XF, and casadi::INTEGRATOR_ZF.

◆ load_integrator()

CASADI_EXPORT void casadi::load_integrator ( const std::string &  name)

Definition at line 126 of file integrator.cpp.

126  {
127  Integrator::load_plugin(name);
128 }

References casadi::PluginInterface< Integrator >::load_plugin().