25 #ifndef CASADI_OPTISTACK_INTERNAL_HPP
26 #define CASADI_OPTISTACK_INTERNAL_HPP
28 #include "optistack.hpp"
29 #include "shared_object.hpp"
42 operator bool()
const {
return ptr_; }
57 friend class InternalOptiCallback;
61 OptiNode(
const std::string& problem_type);
67 MX variable(casadi_int n=1, casadi_int m=1,
const std::string& attribute=
"full");
68 MX variable(
const Sparsity& sp,
const std::string& attribute=
"full");
69 MX variable(
const MX& symbol,
const std::string& attribute=
"full");
72 MX parameter(casadi_int n=1, casadi_int m=1,
const std::string& attribute=
"full");
75 MX parameter(
const Sparsity& sp,
const std::string& attribute=
"full");
77 MX parameter(
const MX& symbol,
const std::string& attribute=
"full");
80 void minimize(
const MX& f,
double linear_scale=1);
83 void subject_to(
const MX& g,
const DM& linear_scale=1,
const Dict& options=
Dict());
88 void solver(
const std::string& solver,
94 void set_initial(
const MX& x,
const DM& v);
95 void set_initial(
const std::vector<MX>& assignments);
104 void set_value(
const MX& x,
const DM& v);
105 void set_value(
const std::vector<MX>& assignments);
109 void set_domain(
const MX& x,
const std::string& domain);
112 void set_linear_scale(
const MX& x,
const DM& scale,
const DM& offset);
115 OptiSol solve(
bool accept_limit);
119 DM value(
const MX& x,
120 const std::vector<MX>& values=std::vector<MX>(),
bool scaled=
false)
const;
122 const std::vector<MX>& values=std::vector<MX>(),
bool scaled=
false)
const {
return x; }
124 const std::vector<MX>& values=std::vector<MX>(),
bool scaled=
false)
const {
136 std::string return_status()
const;
138 bool return_success(
bool accept_limit)
const;
149 std::vector<MX> initial()
const;
152 std::vector<MX> value_variables()
const;
153 std::vector<MX> value_parameters()
const;
156 void callback_class();
157 bool has_callback_class()
const;
160 bool is_parametric(
const MX& expr)
const;
164 std::vector<MX> symvar()
const;
165 std::vector<MX> symvar(
const MX& expr)
const;
170 MetaCon canon_expr(
const MX& expr,
const DM& linear_scale=1)
const;
179 void set_meta(
const MX& m,
const MetaVar& meta);
182 void set_meta_con(
const MX& m,
const MetaCon& meta);
185 void update_user_dict(
const MX& m,
const Dict& meta);
186 Dict user_dict(
const MX& m)
const;
189 MX dual(
const MX& m)
const;
191 void assert_active_symbol(
const MX& m)
const;
196 const std::map<
VariableType, std::vector<DM> >& store)
const;
198 MX x_lookup(casadi_int i)
const;
199 MX g_lookup(casadi_int i)
const;
201 std::string x_describe(casadi_int i,
const Dict& opts=
Dict())
const;
202 std::string g_describe(casadi_int i,
const Dict& opts=
Dict())
const;
203 std::string describe(
const MX& x, casadi_int indent=0,
const Dict& opts=
Dict())
const;
205 void solve_prepare();
206 Function solver_construct(
bool callback=
true);
210 void res(
const DMDict& res);
221 std::string
class_name()
const override {
return "OptiNode"; }
224 casadi_int
nx()
const {
225 if (problem_dirty())
return baked_copy().nx();
226 return nlp_.at(
"x").size1();
230 casadi_int
np()
const {
231 if (problem_dirty())
return baked_copy().np();
232 return nlp_.at(
"p").size1();
236 casadi_int
ng()
const {
237 if (problem_dirty())
return baked_copy().ng();
238 return nlp_.at(
"g").size1();
243 if (problem_dirty())
return baked_copy().x();
249 if (problem_dirty())
return baked_copy().p();
255 if (problem_dirty())
return baked_copy().g();
256 return nlp_unscaled_.at(
"g");
261 if (problem_dirty())
return baked_copy().f();
262 return nlp_unscaled_.at(
"f");
266 if (problem_dirty())
return baked_copy().lbg();
267 return bounds_unscaled_lbg_;
271 if (problem_dirty())
return baked_copy().ubg();
272 return bounds_unscaled_ubg_;
277 if (problem_dirty())
return baked_copy().lam_g();
282 if (problem_dirty())
return baked_copy().x_linear_scale();
283 return DM(linear_scale_);
286 if (problem_dirty())
return baked_copy().x_linear_scale_offset();
287 return DM(linear_scale_offset_);
290 if (problem_dirty())
return baked_copy().g_linear_scale();
291 return DM(g_linear_scale_);
294 if (problem_dirty())
return baked_copy().f_linear_scale();
295 return f_linear_scale_;
297 void assert_empty()
const;
299 void show_infeasibilities(
double tol=0,
const Dict& opts=
Dict())
const;
304 Function to_function(
const std::string& name,
305 const std::vector<MX>& args,
const std::vector<MX>& res,
306 const std::vector<std::string>& name_in,
307 const std::vector<std::string>& name_out,
311 void disp(std::ostream& stream,
bool more=
false)
const override;
316 casadi_int instance_number()
const;
318 static OptiNode* create(
const std::string& problem_type);
332 void assert_solved()
const;
333 void assert_baked()
const;
335 casadi_int g_index_reduce_g(casadi_int i)
const;
336 casadi_int g_index_reduce_x(casadi_int i)
const;
337 casadi_int g_index_unreduce_g(casadi_int i)
const;
341 static std::map<VariableType, std::string> VariableType2String_;
342 std::string variable_type_to_string(
VariableType vt)
const;
344 bool parse_opti_name(
const std::string& name,
VariableType& vt)
const;
345 void register_dual(
MetaCon& meta);
348 void set_value_internal(
const MX& x,
const DM& v,
359 static std::vector<MX> ineq_unchain(
const MX& a,
bool& SWIG_OUTPUT(flipped));
367 const MetaCon& meta_con(
const MX& m)
const;
372 std::vector<MX>
sort(
const std::vector<MX>& v)
const;
375 void assert_has(
const MX& m)
const;
377 bool has(
const MX& m)
const;
380 void assert_has_con(
const MX& m)
const;
382 bool has_con(
const MX& m)
const;
387 std::map<MXNode*, MetaVar> meta_;
389 std::map<MXNode*, MetaCon> meta_con_;
392 std::vector<MX> symbols_;
395 std::vector<bool> discrete_;
400 casadi_int count_var_;
401 casadi_int count_par_;
402 casadi_int count_dual_;
405 std::map< VariableType, std::vector<DM> > store_initial_, store_latest_;
408 std::map< VariableType, std::vector<DM> > store_linear_scale_, store_linear_scale_offset_;
411 std::vector<bool> symbol_active_;
417 mutable std::vector<casadi_int> g_index_reduce_g_;
418 mutable std::vector<casadi_int> g_index_reduce_x_;
419 mutable std::vector<casadi_int> g_index_unreduce_g_;
420 mutable std::vector<casadi_int> target_x_;
421 mutable std::vector<bool> is_simple_;
422 mutable bool reduced_;
430 std::vector<double> linear_scale_;
431 std::vector<double> linear_scale_offset_;
432 std::vector<double> g_linear_scale_;
433 std::vector<double> h_linear_scale_;
435 std::vector<casadi_int> index_all_to_g_;
439 MX bounds_lbg_, bounds_unscaled_lbg_;
440 MX bounds_ubg_, bounds_unscaled_ubg_;
441 std::vector<bool> equality_;
449 double f_linear_scale_;
452 std::string problem_type_;
457 bool old_callback()
const;
459 std::string solver_name_;
460 Dict solver_options_;
462 void assert_only_opti_symbols(
const MX& e)
const;
463 void assert_only_opti_nondual(
const MX& e)
const;
466 static casadi_int instance_count_;
467 casadi_int instance_number_;
470 std::string name_prefix()
const;
472 static std::string format_stacktrace(
const Dict& stacktrace, casadi_int indent);
const Sparsity & sparsity() const
Const access the sparsity - reference to data member.
static Matrix< double > nan(const Sparsity &sp)
create a matrix with all nan
bool problem_dirty() const
void bake()
Fix the structure of the optimization problem.
A simplified interface for NLP modeling/solving.
std::vector< MX > constraints() const
MX g() const
Get all (scalarised) constraint expressions as a column vector.
OptiAdvanced baked_copy() const
DM g_linear_scale() const
MX x() const
Get all (scalarised) decision variables as a symbolic column vector.
casadi_int nx() const
Number of (scalarised) decision variables.
bool solver_dirty() const
casadi_int np() const
Number of (scalarised) parameters.
casadi_int ng() const
Number of (scalarised) constraints.
DM x_linear_scale() const
std::string class_name() const override
Readable name of the internal class.
DM value(const SX &x, const std::vector< MX > &values=std::vector< MX >(), bool scaled=false) const
bool problem_dirty() const
DM x_linear_scale_offset() const
void mark_problem_dirty(bool flag=true)
void mark_solver_dirty(bool flag=true)
MX p() const
Get all (scalarised) parameters as a symbolic column vector.
double f_linear_scale() const
MX lam_g() const
Get dual variables as a symbolic column vector.
MX f() const
Get objective expression.
void mark_solved(bool flag=true)
DM value(const DM &x, const std::vector< MX > &values=std::vector< MX >(), bool scaled=false) const
A simplified interface for NLP modeling/solving.
A simplified interface for NLP modeling/solving.
Pointer that gets set to null when copied.
null_ptr_on_copy(const null_ptr_on_copy &rhs)
std::map< std::string, MX > MXDict
void sort(const std::vector< T > &values, std::vector< T > &sorted_values, std::vector< casadi_int > &indices, bool invert_indices=false)
Sort the data in a vector.
GenericType::Dict Dict
C++ equivalent of Python's dict or MATLAB's struct.
void copy(const T1 *x, casadi_int n, T2 *y)
std::map< std::string, DM > DMDict