Functions
Title

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...
 

Detailed Description

An interpolant function for lookup table data

Parameters
[in]namelabel for the resulting Function
[in]solvername of the plugin
[in]gridcollection of 1D grids whose outer product defines the full N-D rectangular grid
[in]valuesflattened 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])
* 

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
batch_xOT_INTEvaluate a batch of different inputs at once (default 1).casadi::Interpolant
cacheOT_DICTPrepopulate the function cache. Default: emptycasadi::FunctionInternal
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
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
inlineOT_BOOLImplement 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_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
lookup_modeOT_STRINGVECTORSpecifies, 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_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
never_inlineOT_BOOLForbid inlining.casadi::FunctionInternal
output_schemeOT_STRINGVECTORDeprecated option (ignored)casadi::FunctionInternal
post_expandOT_BOOLAfter construction, expand this Function. Default: Falsecasadi::FunctionInternal
post_expand_optionsOT_DICTOptions to be passed to post-construction expansion. Default: emptycasadi::FunctionInternal
print_inOT_BOOLPrint numerical values of inputs [default: false]casadi::FunctionInternal
print_outOT_BOOLPrint numerical values of outputs [default: false]casadi::FunctionInternal
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
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

List of plugins

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


bspline

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


List of available options
IdTypeDescription
algorithmOT_STRINGAlgorithm used for fitting the data: 'not_a_knot' (default, same as Matlab), 'smooth_linear'.
degreeOT_INTVECTORSets, for each grid dimension, the degree of the spline.
linear_solverOT_STRINGSolver used for constructing the coefficient tensor.
linear_solver_optionsOT_DICTOptions to be passed to the linear solver.
smooth_linear_fracOT_DOUBLEWhen 'smooth_linear' algorithm is active, determines sharpness between 0 (sharp, as linear interpolation) and 0.5 (smooth).Default value is 0.1.

linear

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


List of available options
IdTypeDescription
lookup_modeOT_STRINGVECTORSets, 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).
Author
Joel Andersson
Date
2016

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

Function Documentation

◆ doc_interpolant()

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

◆ has_interpolant()

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

◆ interpolant() [1/4]

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

◆ interpolant() [2/4]

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

◆ interpolant() [3/4]

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

◆ interpolant() [4/4]

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

◆ load_interpolant()

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