Functions | |
CASADI_EXPORT Function | casadi::interpolant (const std::string &name, const std::string &solver, const std::vector< std::vector< double > > &grid, const std::vector< double > &values, const Dict &opts=Dict()) |
CASADI_EXPORT Function | casadi::interpolant (const std::string &name, const std::string &solver, const std::vector< std::vector< double > > &grid, casadi_int m=1, const Dict &opts=Dict()) |
Parametric variant of interpolant. More... | |
CASADI_EXPORT Function | casadi::interpolant (const std::string &name, const std::string &solver, const std::vector< casadi_int > &grid_dims, casadi_int m=1, const Dict &opts=Dict()) |
Parametric variant of interpolant. More... | |
CASADI_EXPORT Function | casadi::interpolant (const std::string &name, const std::string &solver, const std::vector< casadi_int > &grid_dims, const std::vector< double > &values, const Dict &opts=Dict()) |
Parametric variant of interpolant. More... | |
CASADI_EXPORT bool | casadi::has_interpolant (const std::string &name) |
Check if a particular plugin is available. More... | |
CASADI_EXPORT void | casadi::load_interpolant (const std::string &name) |
Explicitly load a plugin dynamically. More... | |
CASADI_EXPORT std::string | casadi::doc_interpolant (const std::string &name) |
Get the documentation string for a plugin. More... | |
An interpolant function for lookup table data
[in] | name | label for the resulting Function |
[in] | solver | name of the plugin |
[in] | grid | collection of 1D grids whose outer product defines the full N-D rectangular grid |
[in] | values | flattened vector of all values for all gridpoints |
Syntax 1D
* # Python * xgrid = np.linspace(1,6,6) * V = [-1,-1,-2,-3,0,2] * LUT = casadi.interpolant("LUT","bspline",[xgrid],V) * print(LUT(2.5)) *
* % Matlab * xgrid = 1:6; * V = [-1 -1 -2 -3 0 2]; * LUT = casadi.interpolant('LUT','bspline',{xgrid},V); * LUT(2.5) *
Syntax 2D
* # Python * xgrid = np.linspace(-5,5,11) * ygrid = np.linspace(-4,4,9) * X,Y = np.meshgrid(xgrid,ygrid,indexing='ij') * R = np.sqrt(5*X**2 + Y**2)+ 1 * data = np.sin(R)/R * data_flat = data.ravel(order='F') * LUT = casadi.interpolant('name','bspline',[xgrid,ygrid],data_flat) * print(LUT([0.5,1])) * \enverbatim * \verbatim * % Matlab * xgrid = -5:1:5; * ygrid = -4:1:4; * R = sqrt(5*X.^2 + Y.^2)+ 1; * V = sin(R)./(R); * LUT = interpolant('LUT','bspline',{xgrid, ygrid},V(:)); * LUT([0.5 1]) *
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 |
batch_x | OT_INT | Evaluate a batch of different inputs at once (default 1). | casadi::Interpolant |
cache | OT_DICT | Prepopulate the function cache. Default: empty | casadi::FunctionInternal |
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 |
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 |
inline | OT_BOOL | Implement the lookup table in MX primitives. Useful when you need derivatives with respect to grid and/or coefficients. Such derivatives are fundamentally dense, so use with caution. | casadi::Interpolant |
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 |
lookup_mode | OT_STRINGVECTOR | Specifies, for each grid dimension, the lookup algorithm used to find the correct index. 'linear' uses a for-loop + break; (default when #knots<=100), 'exact' uses floored division (only for uniform grids), 'binary' uses a binary search. (default when #knots>100). | casadi::Interpolant |
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 |
never_inline | OT_BOOL | Forbid inlining. | casadi::FunctionInternal |
output_scheme | OT_STRINGVECTOR | Deprecated option (ignored) | casadi::FunctionInternal |
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_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 |
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 |
- bspline
- linear
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 Interpolant.doc("myextraplugin")
Extra doc: https://github.com/casadi/casadi/wiki/L_239
Id | Type | Description |
---|---|---|
algorithm | OT_STRING | Algorithm used for fitting the data: 'not_a_knot' (default, same as Matlab), 'smooth_linear'. |
degree | OT_INTVECTOR | Sets, for each grid dimension, the degree of the spline. |
linear_solver | OT_STRING | Solver used for constructing the coefficient tensor. |
linear_solver_options | OT_DICT | Options to be passed to the linear solver. |
smooth_linear_frac | OT_DOUBLE | When 'smooth_linear' algorithm is active, determines sharpness between 0 (sharp, as linear interpolation) and 0.5 (smooth).Default value is 0.1. |
Extra doc: https://github.com/casadi/casadi/wiki/L_23b
Id | Type | Description |
---|---|---|
lookup_mode | OT_STRINGVECTOR | Sets, for each grid dimenion, the lookup algorithm used to find the correct index. 'linear' uses a for-loop + break; 'exact' uses floored division (only for uniform grids). |
Extra doc: https://github.com/casadi/casadi/wiki/L_21p
CASADI_EXPORT std::string casadi::doc_interpolant | ( | const std::string & | name | ) |
CASADI_EXPORT bool casadi::has_interpolant | ( | const std::string & | name | ) |
CASADI_EXPORT Function casadi::interpolant | ( | const std::string & | name, |
const std::string & | solver, | ||
const std::vector< casadi_int > & | grid_dims, | ||
casadi_int | m = 1 , |
||
const Dict & | opts = Dict() |
||
) |
The resulting function will have additional arguments for the grid and coefficients
By default, derivatives wrt the coefficients are not supported (zero). Some interpolant plugins may support the inline=true
which enables correct derivatives
Extra doc: https://github.com/casadi/casadi/wiki/L_1p4
CASADI_EXPORT Function casadi::interpolant | ( | const std::string & | name, |
const std::string & | solver, | ||
const std::vector< casadi_int > & | grid_dims, | ||
const std::vector< double > & | values, | ||
const Dict & | opts = Dict() |
||
) |
The resulting function will have an additional argument for the grid
By default, derivatives wrt the coefficients are not supported (zero). Some interpolant plugins may support the inline=true
which enables correct derivatives
Extra doc: https://github.com/casadi/casadi/wiki/L_1p5
CASADI_EXPORT Function casadi::interpolant | ( | const std::string & | name, |
const std::string & | solver, | ||
const std::vector< std::vector< double > > & | grid, | ||
casadi_int | m = 1 , |
||
const Dict & | opts = Dict() |
||
) |
The resulting function will have an additional argument for the coefficients
By default, derivatives wrt the coefficients are not supported (zero). Some interpolant plugins may support the inline=true
which enables correct derivatives
Extra doc: https://github.com/casadi/casadi/wiki/L_1p3
CASADI_EXPORT Function casadi::interpolant | ( | const std::string & | name, |
const std::string & | solver, | ||
const std::vector< std::vector< double > > & | grid, | ||
const std::vector< double > & | values, | ||
const Dict & | opts = Dict() |
||
) |
CASADI_EXPORT void casadi::load_interpolant | ( | const std::string & | name | ) |