Functions
Title

Functions

CASADI_EXPORT Function casadi::conic (const std::string &name, const std::string &solver, const SpDict &qp, const Dict &opts=Dict())
 
CASADI_EXPORT Function casadi::qpsol (const std::string &name, const std::string &solver, const SXDict &qp, const Dict &opts=Dict())
 
CASADI_EXPORT Function casadi::qpsol (const std::string &name, const std::string &solver, const MXDict &qp, const Dict &opts=Dict())
 
CASADI_EXPORT std::vector< std::string > casadi::conic_in ()
 Get input scheme of QP solvers. More...
 
CASADI_EXPORT std::vector< std::string > casadi::conic_out ()
 Get QP solver output scheme of QP solvers. More...
 
CASADI_EXPORT std::string casadi::conic_in (casadi_int ind)
 Get QP solver input scheme name by index. More...
 
CASADI_EXPORT std::string casadi::conic_out (casadi_int ind)
 Get output scheme name by index. More...
 
CASADI_EXPORT casadi_int casadi::conic_n_in ()
 Get the number of QP solver inputs. More...
 
CASADI_EXPORT casadi_int casadi::conic_n_out ()
 Get the number of QP solver outputs. More...
 
CASADI_EXPORT std::vector< std::string > casadi::conic_options (const std::string &name)
 Get all options for a plugin. More...
 
CASADI_EXPORT std::string casadi::conic_option_type (const std::string &name, const std::string &op)
 Get type info for a particular option. More...
 
CASADI_EXPORT std::string casadi::conic_option_info (const std::string &name, const std::string &op)
 Get documentation for a particular option. More...
 
CASADI_EXPORT bool casadi::has_conic (const std::string &name)
 Check if a particular plugin is available. More...
 
CASADI_EXPORT void casadi::load_conic (const std::string &name)
 Explicitly load a plugin dynamically. More...
 
CASADI_EXPORT std::string casadi::doc_conic (const std::string &name)
 Get the documentation string for a plugin. More...
 
CASADI_EXPORT void casadi::conic_debug (const Function &f, const std::string &filename)
 
CASADI_EXPORT void casadi::conic_debug (const Function &f, std::ostream &file)
 

Detailed Description

Create a QP solver Solves the following strictly convex problem:

min          1/2 x' H x + g' x
x

subject to
LBA <= A x <= UBA
LBX <= x   <= UBX

resize(Q x, np, np) + P >= 0 (psd)

with :
H sparse (n x n) positive definite
g dense  (n x 1)
A sparse (nc x n)
Q sparse symmetric (np^2 x n)
P sparse symmetric (np x nq)

n: number of decision variables (x)
nc: number of constraints (A)
nq: shape of psd constraint matrix

If H is not positive-definite, the solver should throw an error.

Second-order cone constraints can be added as psd constraints through a helper function 'soc':

x in R^n y in R

|| x ||_2 <= y

<=>

soc(x, y) psd

This can be proven with soc(x, y)=[y*I x; x' y] using the Shur complement.

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
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
discreteOT_BOOLVECTORIndicates which of the variables are discrete, i.e. integer-valuedcasadi::Conic
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
equalityOT_BOOLVECTORIndicate an upfront hint which of the constraints are equalities. Some solvers may be able to exploit this knowledge. When true, the corresponding lower and upper bounds are assumed equal. When false, the corresponding bounds may be equal or different.casadi::Conic
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
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_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_problemOT_BOOLPrint a numeric description of the problemcasadi::Conic
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

Input scheme: casadi::ConicInput (CONIC_NUM_IN = 12)
Full nameShortDescription
CONIC_HhThe square matrix H: sparse, (n x n). Only the lower triangular part is actually used. The matrix is assumed to be symmetrical.
CONIC_GgThe vector g: dense, (n x 1)
CONIC_AaThe matrix A: sparse, (nc x n) - product with x must be dense.
CONIC_LBAlbadense, (nc x 1)
CONIC_UBAubadense, (nc x 1)
CONIC_LBXlbxdense, (n x 1)
CONIC_UBXubxdense, (n x 1)
CONIC_X0x0dense, (n x 1)
CONIC_LAM_X0lam_x0dense
CONIC_LAM_A0lam_a0dense
CONIC_QqThe matrix Q: sparse symmetric, (np^2 x n)
CONIC_PpThe matrix P: sparse symmetric, (np x np)

Output scheme: casadi::ConicOutput (CONIC_NUM_OUT = 4)
Full nameShortDescription
CONIC_XxThe primal solution.
CONIC_COSTcostThe optimal cost.
CONIC_LAM_Alam_aThe dual solution corresponding to linear bounds.
CONIC_LAM_Xlam_xThe dual solution corresponding to simple bounds.

List of plugins

- cbc

- clp

- cplex

- daqp

- fatrop

- gurobi

- highs

- hpipm

- hpmpc

- ooqp

- osqp

- proxqp

- qpoases

- sqic

- superscs

- ipqp

- nlpsol

- qrqp

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 Conic.doc("myextraplugin")


cbc

Interface to Cbc solver for sparse Quadratic Programs

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


List of available options
IdTypeDescription
cbcOT_DICTOptions to be passed to CBC.Three sets of options are supported. The first can be found in OsiSolverParameters.hpp. The second can be found in CbcModel.hpp. The third are options that can be passed to CbcMain1.
hot_startOT_BOOLHot start with x0 [Default false].
sos_groupsOT_INTVECTORVECTORDefinition of SOS groups by indices.
sos_typesOT_INTVECTORSpecify 1 or 2 for each SOS group.
sos_weightsOT_DOUBLEVECTORVECTORWeights corresponding to SOS entries.

clp

Interface to Clp solver for sparse Quadratic Programs

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


List of available options
IdTypeDescription
clpOT_DICTOptions to be passed to CLP. A first set of options can be found in ClpParameters.hpp. eg. 'PrimalTolerance'. There are other options in additions. 'AutomaticScaling' (bool) is recognised. 'initial_solve' (default off) activates the use of Clp's initialSolve. 'initial_solve_options' takes a dictionary with following keys (see ClpSolve.hpp): SolveType (string), PresolveType (string), NumberPasses, SpecialOptions (intvectorvector), IndependentOptions (intvectorvector).

cplex

Interface to Cplex solver for sparse Quadratic Programs

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


List of available options
IdTypeDescription
cplexOT_DICTOptions to be passed to CPLEX
dep_checkOT_INTDetect redundant constraints.
dump_filenameOT_STRINGThe filename to dump to.
dump_to_fileOT_BOOLDumps QP to file in CPLEX format.
mip_startOT_BOOLHot start integers with x0 [Default false].
qp_methodOT_INTDetermines which CPLEX algorithm to use.
sos_groupsOT_INTVECTORVECTORDefinition of SOS groups by indices.
sos_typesOT_INTVECTORSpecify 1 or 2 for each SOS group.
sos_weightsOT_DOUBLEVECTORVECTORWeights corresponding to SOS entries.
tolOT_DOUBLETolerance of solver
version_suffixOT_STRINGSpecify version of cplex to load. We will attempt to load libcplex<version_suffix>.[so|dll|dylib]. Default value is taken from CPLEX_VERSION env variable.
warm_startOT_BOOLUse warm start with simplex methods (affects only the simplex methods).

daqp

Interface to Daqp solver for sparse Quadratic Programs, see daqp.dev for more information and https://www.maths.ed.ac.uk/hall/Daqp/DaqpOptions.html for a list of options.

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


List of available options
IdTypeDescription
daqpOT_DICTOptions to be passed to Daqp.

fatrop

Interface to Fatrop Solver

In order to use this interface, you must:

       A0 B0 -I
       C0 D0
              A1 B1 -I
              C1 D1

where I must be a diagonal sparse matrix


List of available options
IdTypeDescription
NOT_INTOCP horizon
fatropOT_DICTOptions to be passed to fatrop
ngOT_INTVECTORNumber of non-dynamic constraints, length N+1
nuOT_INTVECTORNumber of controls, length N
nxOT_INTVECTORNumber of states, length N+1
structure_detectionOT_STRINGNONE | auto | manual

gurobi

Interface to the GUROBI Solver for quadratic programming

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


List of available options
IdTypeDescription
gurobiOT_DICTOptions to be passed to gurobi.
sos_groupsOT_INTVECTORVECTORDefinition of SOS groups by indices.
sos_typesOT_INTVECTORSpecify 1 or 2 for each SOS group.
sos_weightsOT_DOUBLEVECTORVECTORWeights corresponding to SOS entries.
vtypeOT_STRINGVECTORType of variables: [CONTINUOUS|binary|integer|semicont|semiint]

highs

Interface to HiGHS solver for sparse Quadratic Programs, see highs.dev for more information and https://www.maths.ed.ac.uk/hall/HiGHS/HighsOptions.html for a list of options.

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


List of available options
IdTypeDescription
highsOT_DICTOptions to be passed to HiGHS.

hpipm

Interface to HPIPM Solver

In order to use this interface, you must:

       A0 B0 -I
       C0 D0
              A1 B1 -I
              C1 D1

where I must be a diagonal sparse matrix


List of available options
IdTypeDescription
NOT_INTOCP horizon
hpipmOT_DICTOptions to be passed to hpipm
infOT_DOUBLEReplace infinities by this amount [default: 1e8]
ngOT_INTVECTORNumber of non-dynamic constraints, length N+1
nuOT_INTVECTORNumber of controls, length N
nxOT_INTVECTORNumber of states, length N+1

hpmpc

Interface to HMPC Solver

In order to use this interface, you must:

       A0 B0 -I
       C0 D0
              A1 B1 -I
              C1 D1

where I must be a diagonal sparse matrix


List of available options
IdTypeDescription
NOT_INTOCP horizon
blasfeo_targetOT_STRINGhpmpc target
infOT_DOUBLEHPMPC cannot handle infinities. Infinities will be replaced by this option's value.
max_iterOT_INTMax number of iterations
mu0OT_DOUBLEMax element in cost function as estimate of max multiplier
ngOT_INTVECTORNumber of non-dynamic constraints, length N+1
nuOT_INTVECTORNumber of controls, length N
nxOT_INTVECTORNumber of states, length N+1
print_levelOT_INTAmount of diagnostic printing [Default: 1].
targetOT_STRINGhpmpc target
tolOT_DOUBLETolerance in the duality measure
warm_startOT_BOOLUse warm-starting

ooqp

Interface to the OOQP Solver for quadratic programming The current implementation assumes that OOQP is configured with the MA27 sparse linear solver.

NOTE: when doing multiple calls to evaluate(), check if you need to reInit();

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


List of available options
IdTypeDescription
artolOT_DOUBLEtolerance as provided with setArTol to OOQP
mutolOT_DOUBLEtolerance as provided with setMuTol to OOQP
print_levelOT_INTPrint level. OOQP listens to print_level 0, 10 and 100

osqp

Interface to the OSQP Solver for quadratic programming

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

Interface to the PROXQP Solver for quadratic programming

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


List of available options
IdTypeDescription
osqpOT_DICTconst Options to be passed to osqp.
warm_start_dualOT_BOOLUse lam_a0 and lam_x0 input to warmstart [Default: truw].
warm_start_primalOT_BOOLUse x0 input to warmstart [Default: true].

proxqp


List of available options
IdTypeDescription
proxqpOT_DICTconst proxqp options.
warm_start_dualOT_BOOLUse y and z input to warmstart [Default: true].
warm_start_primalOT_BOOLUse x input to warmstart [Default: true].

qpoases

Interface to QPOases Solver for quadratic programming

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


List of available options
IdTypeDescription
CPUtimeOT_DOUBLEThe maximum allowed CPU time in seconds for the whole initialisation (and the actually required one on output). Disabled if unset.
boundRelaxationOT_DOUBLEInitial relaxation of bounds to start homotopy and initial value for far bounds.
boundToleranceOT_DOUBLEIf upper and lower bounds differ less than this tolerance, they are regarded equal, i.e. as equality constraint.
enableCholeskyRefactorisationOT_INTSpecifies the frequency of a full re-factorisation of projected Hessian matrix: 0: turns them off, 1: uses them at each iteration etc.
enableDriftCorrectionOT_INTSpecifies the frequency of drift corrections: 0: turns them off.
enableEqualitiesOT_BOOLSpecifies whether equalities should be treated as always active (True) or not (False)
enableFarBoundsOT_BOOLEnables the use of far bounds.
enableFlippingBoundsOT_BOOLEnables the use of flipping bounds.
enableFullLITestsOT_BOOLEnables condition-hardened (but more expensive) LI test.
enableInertiaCorrectionOT_BOOLShould working set be repaired when negative curvature is discovered during hotstart.
enableNZCTestsOT_BOOLEnables nonzero curvature tests.
enableRampingOT_BOOLEnables ramping.
enableRegularisationOT_BOOLEnables automatic Hessian regularisation.
epsDenOT_DOUBLEDenominator tolerance for ratio tests.
epsFlippingOT_DOUBLETolerance of squared Cholesky diagonal factor which triggers flipping bound.
epsIterRefOT_DOUBLEEarly termination tolerance for iterative refinement.
epsLITestsOT_DOUBLETolerance for linear independence tests.
epsNZCTestsOT_DOUBLETolerance for nonzero curvature tests.
epsNumOT_DOUBLENumerator tolerance for ratio tests.
epsRegularisationOT_DOUBLEScaling factor of identity matrix used for Hessian regularisation.
finalRampingOT_DOUBLEFinal value for ramping strategy.
growFarBoundsOT_DOUBLEFactor to grow far bounds.
hessian_typeOT_STRINGType of Hessian - see qpOASES documentation [UNKNOWN|posdef|semidef|indef|zero|identity]]
initialFarBoundsOT_DOUBLEInitial size for far bounds.
initialRampingOT_DOUBLEStart value for ramping strategy.
initialStatusBoundsOT_STRINGInitial status of bounds at first iteration.
linsol_pluginOT_STRINGLinear solver plugin
maxDualJumpOT_DOUBLEMaximum allowed jump in dual variables in linear independence tests.
maxPrimalJumpOT_DOUBLEMaximum allowed jump in primal variables in nonzero curvature tests.
max_schurOT_INTMaximal number of Schur updates [75]
nWSROT_INTThe maximum number of working set recalculations to be performed during the initial homotopy. Default is 5(nx + nc)
numRefinementStepsOT_INTMaximum number of iterative refinement steps.
numRegularisationStepsOT_INTMaximum number of successive regularisation steps.
printLevelOT_STRINGDefines the amount of text output during QP solution, see Section 5.7
schurOT_BOOLUse Schur Complement Approach [false]
sparseOT_BOOLFormulate the QP using sparse matrices. [false]
terminationToleranceOT_DOUBLERelative termination tolerance to stop homotopy.

sqic

Interface to the SQIC solver for quadratic programming

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


superscs

Interface to the SuperSCS solver for conic programming

Joris Gillis, 2019

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


List of available options
IdTypeDescription
superscsOT_DICTOptions to be passed to superscs.

ipqp

Solves QPs using a Mehrotra predictor-corrector interior point method

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


List of available options
IdTypeDescription
constr_viol_tolOT_DOUBLEConstraint violation tolerance [1e-8].
dual_inf_tolOT_DOUBLEDual feasibility violation tolerance [1e-8]
linear_solverOT_STRINGA custom linear solver creator function [default: ldl]
linear_solver_optionsOT_DICTOptions to be passed to the linear solver
max_iterOT_INTMaximum number of iterations [1000].
min_lamOT_DOUBLESmallest multiplier treated as inactive for the initial active set [0].
print_headerOT_BOOLPrint header [true].
print_infoOT_BOOLPrint info [true].
print_iterOT_BOOLPrint iterations [true].

nlpsol

Solve QPs using an Nlpsol Use the 'nlpsol' option to specify the NLP solver to use.

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


List of available options
IdTypeDescription
nlpsolOT_STRINGName of solver.
nlpsol_optionsOT_DICTOptions to be passed to solver.

qrqp

Solve QPs using an active-set method

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


List of available options
IdTypeDescription
constr_viol_tolOT_DOUBLEConstraint violation tolerance [1e-8].
dual_inf_tolOT_DOUBLEDual feasibility violation tolerance [1e-8]
max_iterOT_INTMaximum number of iterations [1000].
min_lamOT_DOUBLESmallest multiplier treated as inactive for the initial active set [0].
print_headerOT_BOOLPrint header [true].
print_infoOT_BOOLPrint info [true].
print_iterOT_BOOLPrint iterations [true].
print_lincombOT_BOOLPrint dependant linear combinations of constraints [false]. Printed numbers are 0-based indices into the vector of [simple bounds;linear bounds]
Author
Joel Andersson
Date
2011-2015

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

Function Documentation

◆ conic()

CASADI_EXPORT Function casadi::conic ( const std::string &  name,
const std::string &  solver,
const SpDict qp,
const Dict opts = Dict() 
)

◆ conic_debug() [1/2]

CASADI_EXPORT void casadi::conic_debug ( const Function f,
const std::string &  filename 
)

Generate native code in the interfaced language for debugging

◆ conic_debug() [2/2]

CASADI_EXPORT void casadi::conic_debug ( const Function f,
std::ostream &  file 
)

Generate native code in the interfaced language for debugging

◆ conic_in() [1/2]

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

◆ conic_in() [2/2]

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

◆ conic_n_in()

CASADI_EXPORT casadi_int casadi::conic_n_in ( )

◆ conic_n_out()

CASADI_EXPORT casadi_int casadi::conic_n_out ( )

◆ conic_option_info()

CASADI_EXPORT std::string casadi::conic_option_info ( const std::string &  name,
const std::string &  op 
)

◆ conic_option_type()

CASADI_EXPORT std::string casadi::conic_option_type ( const std::string &  name,
const std::string &  op 
)

◆ conic_options()

CASADI_EXPORT std::vector<std::string> casadi::conic_options ( const std::string &  name)

◆ conic_out() [1/2]

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

◆ conic_out() [2/2]

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

◆ doc_conic()

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

◆ has_conic()

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

◆ load_conic()

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

◆ qpsol() [1/2]

CASADI_EXPORT Function casadi::qpsol ( const std::string &  name,
const std::string &  solver,
const MXDict qp,
const Dict opts = Dict() 
)

◆ qpsol() [2/2]

CASADI_EXPORT Function casadi::qpsol ( const std::string &  name,
const std::string &  solver,
const SXDict qp,
const Dict opts = Dict() 
)