Functions
Title

Functions

Function casadi::dplesol (const std::string &name, const std::string &solver, const SpDict &st, const Dict &opts)
 
MX casadi::dplesol (const MX &A, const MX &V, const std::string &solver, const Dict &opts)
 
CASADI_EXPORT MXVector casadi::dplesol (const MXVector &A, const MXVector &V, const std::string &solver, const Dict &opts)
 
CASADI_EXPORT DMVector casadi::dplesol (const DMVector &A, const DMVector &V, const std::string &solver, const Dict &opts)
 
std::vector< std::string > casadi::dple_in ()
 Get input scheme of DPLE solvers. More...
 
std::vector< std::string > casadi::dple_out ()
 Get output scheme of DPLE solvers. More...
 
std::string casadi::dple_in (casadi_int ind)
 Get DPLE input scheme name by index. More...
 
std::string casadi::dple_out (casadi_int ind)
 Get DPLE output scheme name by index. More...
 
casadi_int casadi::dple_n_in ()
 Get the number of QP solver inputs. More...
 
casadi_int casadi::dple_n_out ()
 Get the number of QP solver outputs. More...
 
bool casadi::has_dple (const std::string &name)
 Check if a particular plugin is available. More...
 
void casadi::load_dple (const std::string &name)
 Explicitly load a plugin dynamically. More...
 
std::string casadi::doc_dple (const std::string &name)
 Get the documentation string for a plugin. More...
 

Detailed Description

Discrete periodic Lyapunov Equation solver Given matrices $A_k$ and symmetric $V_k, k = 0..K-1$

A_k in R^(n x n)
V_k in R^n

provides all of $P_k$ that satisfy:

P_0 = A_(K-1)*P_(K-1)*A_(K-1)' + V_k
P_k+1 = A_k*P_k*A_k' + V_k  for k = 1..K-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
cacheOT_DICTPrepopulate the function cache. Default: emptycasadi::FunctionInternal
compilerOT_STRINGJust-in-time compiler plugin to be used.casadi::FunctionInternal
const_dimOT_BOOLAssume constant dimension of Pcasadi::Dple
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
eps_unstableOT_DOUBLEA margin for unstability detectioncasadi::Dple
error_on_failOT_BOOLThrow exceptions when function evaluation fails (default true).casadi::ProtoFunction
error_unstableOT_BOOLThrow an exception when it is detected that Product(A_i, i=N..1)has eigenvalues greater than 1-eps_unstablecasadi::Dple
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
pos_defOT_BOOLAssume P positive definitecasadi::Dple
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

Input scheme: casadi::DpleInput (DPLE_NUM_IN = 2)
Full nameShortDescription
DPLE_AaA matrices (horzcat when const_dim, diagcat otherwise) [a].
DPLE_VvV matrices (horzcat when const_dim, diagcat otherwise) [v].

Output scheme: casadi::DpleOutput (DPLE_NUM_OUT = 1)
Full nameShortDescription
DPLE_PpLyapunov matrix (horzcat when const_dim, diagcat otherwise) (Cholesky of P if pos_def) [p].

List of plugins

- slicot

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


slicot

An efficient solver for Discrete Periodic Lyapunov Equations using SLICOT

Uses Periodic Schur Decomposition ('psd') and does not assume positive definiteness. Based on Periodic Lyapunov equations: some applications and new algorithms. Int. J. Control, vol. 67, pp. 69-87, 1997.

Overview of the method: J. Gillis Practical Methods for Approximate Robust Periodic Optimal Control ofNonlinear Mechanical Systems, PhD Thesis, KULeuven, 2015

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


List of available options
IdTypeDescription
linear_solverOT_STRINGUser-defined linear solver class. Needed for sensitivities.
linear_solver_optionsOT_DICTOptions to be passed to the linear solver.
psd_num_zeroOT_DOUBLENumerical zero used in Periodic Schur decomposition with slicot.This option is needed when your systems has Floquet multiplierszero or close to zero
Author
Joris Gillis
Date
2013-2016

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

Function Documentation

◆ doc_dple()

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

Definition at line 39 of file dple.cpp.

39  {
40  return Dple::getPlugin(name).doc;
41  }

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

◆ dple_in() [1/2]

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

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

Definition at line 102 of file dple.cpp.

102  {
103  std::vector<std::string> ret(dple_n_in());
104  for (size_t i=0; i<ret.size(); ++i) ret[i]=dple_in(i);
105  return ret;
106  }
casadi_int dple_n_in()
Get the number of QP solver inputs.
Definition: dple.cpp:131
std::string dple_in(casadi_int ind)
Get DPLE input scheme name by index.
Definition: dple.cpp:114

References casadi::dple_n_in().

Referenced by casadi::Dple::get_name_in().

◆ dple_in() [2/2]

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

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

Definition at line 114 of file dple.cpp.

114  {
115  switch (static_cast<DpleInput>(ind)) {
116  case DPLE_A: return "a";
117  case DPLE_V: return "v";
118  case DPLE_NUM_IN: break;
119  }
120  return std::string();
121  }
DpleInput
Input arguments of a dple solver [dpleIn].
Definition: dple.hpp:119
@ DPLE_NUM_IN
Definition: dple.hpp:124
@ DPLE_A
A matrices (horzcat when const_dim, diagcat otherwise) [a].
Definition: dple.hpp:121
@ DPLE_V
V matrices (horzcat when const_dim, diagcat otherwise) [v].
Definition: dple.hpp:123

References casadi::DPLE_A, casadi::DPLE_NUM_IN, and casadi::DPLE_V.

◆ dple_n_in()

CASADI_EXPORT casadi_int casadi::dple_n_in ( )

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

Definition at line 131 of file dple.cpp.

131  {
132  return DPLE_NUM_IN;
133  }

References casadi::DPLE_NUM_IN.

Referenced by casadi::dple_in().

◆ dple_n_out()

CASADI_EXPORT casadi_int casadi::dple_n_out ( )

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

Definition at line 135 of file dple.cpp.

135  {
136  return DPLE_NUM_OUT;
137  }
@ DPLE_NUM_OUT
Number of arguments.
Definition: dple.hpp:132

References casadi::DPLE_NUM_OUT.

Referenced by casadi::dple_out().

◆ dple_out() [1/2]

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

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

Definition at line 108 of file dple.cpp.

108  {
109  std::vector<std::string> ret(dple_n_out());
110  for (size_t i=0; i<ret.size(); ++i) ret[i]=dple_out(i);
111  return ret;
112  }
std::string dple_out(casadi_int ind)
Get DPLE output scheme name by index.
Definition: dple.cpp:123
casadi_int dple_n_out()
Get the number of QP solver outputs.
Definition: dple.cpp:135

References casadi::dple_n_out().

Referenced by casadi::Dple::get_name_out().

◆ dple_out() [2/2]

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

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

Definition at line 123 of file dple.cpp.

123  {
124  switch (static_cast<DpleOutput>(ind)) {
125  case DPLE_P: return "p";
126  case DPLE_NUM_OUT: break;
127  }
128  return std::string();
129  }
DpleOutput
Output arguments of a dple solver [dpleOut].
Definition: dple.hpp:128
@ DPLE_P
Lyapunov matrix (horzcat when const_dim, diagcat otherwise) (Cholesky of P if pos_def) [p].
Definition: dple.hpp:130

References casadi::DPLE_NUM_OUT, and casadi::DPLE_P.

◆ dplesol() [1/4]

CASADI_EXPORT DMVector casadi::dplesol ( const DMVector A,
const DMVector V,
const std::string &  solver,
const Dict opts 
)

Definition at line 71 of file dple.cpp.

72  {
73  casadi_assert(A.size()==V.size(),
74  "dplesol: sizes of A vector (" + str(A.size()) + ") and V vector "
75  "(" + str(V.size()) + ") must match.");
76  std::vector<DM> Adense, Vdense;
77 
78  for (casadi_int i=0;i<A.size();++i) {
79  Adense.push_back(densify(A[i]));
80  Vdense.push_back(densify(V[i]));
81  }
82 
83  DM Afull = diagcat(Adense);
84  DM Vfull = diagcat(Vdense);
85 
86  SpDict sp;
87  sp["a"] = Afull.sparsity();
88  sp["v"] = Vfull.sparsity();
89  Function f = dplesol("dplesol", solver, sp, opts);
90  DMDict f_in;
91  f_in["a"] = Afull;
92  f_in["v"] = Vfull;
93  DMDict f_out = f(f_in);
94  return diagsplit(f_out["p"], f_out["p"].size1()/A.size());
95  }
Function dplesol(const std::string &name, const std::string &solver, const SpDict &st, const Dict &opts)
Definition: dple.cpp:97
std::string str(const T &v)
String representation, any type.
std::map< std::string, Sparsity > SpDict
Definition: sparsity.hpp:1502
Matrix< double > DM
Definition: dm_fwd.hpp:33
std::map< std::string, DM > DMDict
Definition: dm_fwd.hpp:36

References casadi::dplesol(), casadi::Matrix< Scalar >::sparsity(), and casadi::str().

◆ dplesol() [2/4]

CASADI_EXPORT MX casadi::dplesol ( const MX A,
const MX V,
const std::string &  solver,
const Dict opts 
)

Definition at line 43 of file dple.cpp.

43  {
44  SpDict sp;
45  sp["a"] = A.sparsity();
46  sp["v"] = V.sparsity();
47  Function f = dplesol("dplesol", solver, sp, opts);
48  MXDict f_in;
49  f_in["a"] = A;
50  f_in["v"] = V;
51  MXDict f_out = f(f_in);
52  return f_out["p"];
53  }
std::map< std::string, MX > MXDict
Definition: mx.hpp:1009

References casadi::MX::sparsity().

Referenced by casadi::dplesol(), casadi::Dple::get_forward(), and casadi::Dple::get_reverse().

◆ dplesol() [3/4]

CASADI_EXPORT MXVector casadi::dplesol ( const MXVector A,
const MXVector V,
const std::string &  solver,
const Dict opts 
)

Definition at line 55 of file dple.cpp.

56  {
57  casadi_assert(A.size()==V.size(),
58  "dplesol: sizes of A vector (" + str(A.size()) + ") and V vector "
59  "(" + str(V.size()) + ") must match.");
60  std::vector<MX> Adense, Vdense;
61 
62  for (casadi_int i=0;i<A.size();++i) {
63  Adense.push_back(densify(A[i]));
64  Vdense.push_back(densify(V[i]));
65  }
66 
67  MX ret = dplesol(diagcat(Adense), diagcat(Vdense), solver, opts);
68  return diagsplit(ret, ret.size1()/A.size());
69  }

References casadi::dplesol(), casadi::GenericMatrix< MatType >::size1(), and casadi::str().

◆ dplesol() [4/4]

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

Definition at line 97 of file dple.cpp.

98  {
99  return Function::create(Dple::instantiate(name, solver, st), opts);
100  }

References casadi::Function::create(), and casadi::PluginInterface< Dple >::instantiate().

◆ has_dple()

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

Definition at line 31 of file dple.cpp.

31  {
32  return Dple::has_plugin(name);
33  }

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

◆ load_dple()

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

Definition at line 35 of file dple.cpp.

35  {
36  Dple::load_plugin(name);
37  }

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