optistack_internal.hpp
1 /*
2  * This file is part of CasADi.
3  *
4  * CasADi -- A symbolic framework for dynamic optimization.
5  * Copyright (C) 2010-2023 Joel Andersson, Joris Gillis, Moritz Diehl,
6  * KU Leuven. All rights reserved.
7  * Copyright (C) 2011-2014 Greg Horn
8  *
9  * CasADi is free software; you can redistribute it and/or
10  * modify it under the terms of the GNU Lesser General Public
11  * License as published by the Free Software Foundation; either
12  * version 3 of the License, or (at your option) any later version.
13  *
14  * CasADi is distributed in the hope that it will be useful,
15  * but WITHOUT ANY WARRANTY; without even the implied warranty of
16  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
17  * Lesser General Public License for more details.
18  *
19  * You should have received a copy of the GNU Lesser General Public
20  * License along with CasADi; if not, write to the Free Software
21  * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA
22  *
23  */
24 
25 #ifndef CASADI_OPTISTACK_INTERNAL_HPP
26 #define CASADI_OPTISTACK_INTERNAL_HPP
27 
28 #include "optistack.hpp"
29 #include "shared_object_internal.hpp"
30 
31 namespace casadi {
32 
33 #ifndef SWIG
35 template<class T>
36 class null_ptr_on_copy {
37 public:
38  null_ptr_on_copy<T>() : ptr_(nullptr) {}
39  null_ptr_on_copy<T>(const null_ptr_on_copy<T>& rhs) : ptr_(nullptr) {}
40  void operator=(T* ptr) { ptr_ = ptr; }
41  T* operator->() { return ptr_; }
42  operator bool() const { return ptr_; }
43 private:
44  T* ptr_;
45 };
46 #endif
47 
55 class CASADI_EXPORT OptiNode :
56  public SharedObjectInternal {
57  friend class InternalOptiCallback;
58 public:
59 
61  OptiNode(const std::string& problem_type);
62 
65 
67  MX variable(casadi_int n=1, casadi_int m=1, const std::string& attribute="full");
68 
70  MX parameter(casadi_int n=1, casadi_int m=1, const std::string& attribute="full");
71 
73  void minimize(const MX& f);
74 
76  void subject_to(const MX& g);
78  void subject_to();
79 
81  void solver(const std::string& solver,
82  const Dict& plugin_options=Dict(),
83  const Dict& solver_options=Dict());
84 
87  void set_initial(const MX& x, const DM& v);
88  void set_initial(const std::vector<MX>& assignments);
90 
92 
97  void set_value(const MX& x, const DM& v);
98  void set_value(const std::vector<MX>& assignments);
100 
102  void set_domain(const MX& x, const std::string& domain);
103 
105  OptiSol solve(bool accept_limit);
106 
109  DM value(const MX& x, const std::vector<MX>& values=std::vector<MX>()) const;
110  DM value(const DM& x, const std::vector<MX>& values=std::vector<MX>()) const { return x; }
111  DM value(const SX& x, const std::vector<MX>& values=std::vector<MX>()) const {
112  return DM::nan(x.sparsity());
113  }
115 
117  Opti copy() const;
118 
120  Dict stats() const;
121 
123  std::string return_status() const;
125  bool return_success(bool accept_limit) const;
126 
129 
131  std::vector<MX> initial() const;
132 
134  std::vector<MX> value_variables() const;
135  std::vector<MX> value_parameters() const;
136 
137  void callback_class(OptiCallback* callback);
139  bool has_callback_class() const;
140 
142  bool is_parametric(const MX& expr) const;
143 
146  std::vector<MX> symvar() const;
147  std::vector<MX> symvar(const MX& expr) const;
148  std::vector<MX> symvar(const MX& expr, VariableType type) const;
150 
152  MetaCon canon_expr(const MX& expr) const;
153 
155  MetaVar get_meta(const MX& m) const;
156 
158  MetaCon get_meta_con(const MX& m) const;
159 
161  void set_meta(const MX& m, const MetaVar& meta);
162 
164  void set_meta_con(const MX& m, const MetaCon& meta);
165 
167  void update_user_dict(const MX& m, const Dict& meta);
168  Dict user_dict(const MX& m) const;
169 
171  MX dual(const MX& m) const;
172 
173  void assert_active_symbol(const MX& m) const;
174 
175  std::vector<MX> active_symvar(VariableType type) const;
176  std::vector<DM> active_values(VariableType type) const;
177 
178  MX x_lookup(casadi_int i) const;
179  MX g_lookup(casadi_int i) const;
180 
181  std::string x_describe(casadi_int i) const;
182  std::string g_describe(casadi_int i) const;
183  std::string describe(const MX& x, casadi_int indent=0) const;
184 
186  Function solver_construct(bool callback=true);
187  DMDict solve_actual(const DMDict& args);
188 
189  DMDict arg() const { return arg_; }
190  void res(const DMDict& res);
191  DMDict res() const { return res_; }
192  std::vector<MX> constraints() const { return g_; }
193  MX objective() const { return f_; }
194 
196  OptiAdvanced s = copy();
197  if (s.problem_dirty()) s.bake();
198  return s;
199  }
200 
201  std::string class_name() const override { return "OptiNode"; }
202 
204  casadi_int nx() const {
205  if (problem_dirty()) return baked_copy().nx();
206  return nlp_.at("x").size1();
207  }
208 
210  casadi_int np() const {
211  if (problem_dirty()) return baked_copy().np();
212  return nlp_.at("p").size1();
213  }
214 
216  casadi_int ng() const {
217  if (problem_dirty()) return baked_copy().ng();
218  return nlp_.at("g").size1();
219  }
220 
222  MX x() const {
223  if (problem_dirty()) return baked_copy().x();
224  return nlp_.at("x");
225  }
226 
228  MX p() const {
229  if (problem_dirty()) return baked_copy().p();
230  return nlp_.at("p");
231  }
232 
234  MX g() const {
235  if (problem_dirty()) return baked_copy().g();
236  return nlp_.at("g");
237  }
238 
240  MX f() const {
241  if (problem_dirty()) return baked_copy().f();
242  return nlp_.at("f");
243  }
244 
245  MX lbg() const {
246  if (problem_dirty()) return baked_copy().lbg();
247  return bounds_lbg_;
248  }
249 
250  MX ubg() const {
251  if (problem_dirty()) return baked_copy().ubg();
252  return bounds_ubg_;
253  }
254 
256  MX lam_g() const {
257  if (problem_dirty()) return baked_copy().lam_g();
258  return lam_;
259  }
260  void assert_empty() const;
261 
262 
266  Function to_function(const std::string& name,
267  const std::vector<MX>& args, const std::vector<MX>& res,
268  const std::vector<std::string>& name_in,
269  const std::vector<std::string>& name_out,
270  const Dict& opts);
271 
273  void disp(std::ostream& stream, bool more=false) const override;
274 
276  void bake();
277 
278  casadi_int instance_number() const;
279 
280  static OptiNode* create(const std::string& problem_type);
281 
283  void mark_problem_dirty(bool flag=true) { problem_dirty_=flag; mark_solver_dirty(); }
284  bool problem_dirty() const { return problem_dirty_; }
285 
287  void mark_solver_dirty(bool flag=true) { solver_dirty_=flag; mark_solved(false); }
288  bool solver_dirty() const { return solver_dirty_; }
289 
290  bool solved_;
291  void mark_solved(bool flag=true) { solved_ = flag;}
292  bool solved() const { return solved_; }
293 
294  void assert_solved() const;
295  void assert_baked() const;
296 
297 private:
298 
299  static std::map<VariableType, std::string> VariableType2String_;
300  std::string variable_type_to_string(VariableType vt) const;
301 
302  bool parse_opti_name(const std::string& name, VariableType& vt) const;
303  void register_dual(MetaCon& meta);
304 
306  void set_value_internal(const MX& x, const DM& v);
307 
316  static std::vector<MX> ineq_unchain(const MX& a, bool& SWIG_OUTPUT(flipped));
317 
319  const MetaVar& meta(const MX& m) const;
321  MetaVar& meta(const MX& m);
322 
324  const MetaCon& meta_con(const MX& m) const;
326  MetaCon& meta_con(const MX& m);
327 
329  std::vector<MX> sort(const std::vector<MX>& v) const;
330 
332  void assert_has(const MX& m) const;
333 
334  bool has(const MX& m) const;
335 
337  void assert_has_con(const MX& m) const;
338 
339  bool has_con(const MX& m) const;
340 
341  // Data members
342 
344  std::map<MXNode*, MetaVar> meta_;
346  std::map<MXNode*, MetaCon> meta_con_;
347 
349  std::vector<MX> symbols_;
350 
351  // Which x entries are discrete?
352  std::vector<bool> discrete_;
353 
355  casadi_int count_;
356 
357  casadi_int count_var_;
358  casadi_int count_par_;
359  casadi_int count_dual_;
360 
362  std::map< VariableType, std::vector<DM> > store_initial_, store_latest_;
363 
365  std::vector<bool> symbol_active_;
366 
368  Function solver_;
369 
371  DMDict res_;
372  DMDict arg_;
373  MXDict nlp_;
374  MX lam_;
375 
377  Function bounds_;
378  MX bounds_lbg_;
379  MX bounds_ubg_;
380  std::vector<bool> equality_;
381 
383  std::vector<MX> g_;
384 
386  MX f_;
387 
389  std::string problem_type_;
390 
391  null_ptr_on_copy<OptiCallback> user_callback_;
392  Function callback_;
393 
394  bool old_callback() const;
395 
396  std::string solver_name_;
397  Dict solver_options_;
398 
399  void assert_only_opti_symbols(const MX& e) const;
400  void assert_only_opti_nondual(const MX& e) const;
401 
402 
403  static casadi_int instance_count_;
404  casadi_int instance_number_;
405 
406 
407  std::string name_prefix() const;
408 
409  static std::string format_stacktrace(const Dict& stacktrace, casadi_int indent);
410 
411 };
412 
413 
414 } // namespace casadi
415 
416 #endif // CASADI_OPTISTACK_INTERNAL_HPP
Function object.
Definition: function.hpp:60
Sparsity sparsity() const
Get the sparsity pattern.
MX - Matrix expression.
Definition: mx.hpp:84
Sparse matrix class. SX and DM are specializations.
Definition: matrix_decl.hpp:92
static Matrix< Scalar > 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.
Function solver_construct(bool callback=true)
std::string g_describe(casadi_int i) const
Dict user_dict(const MX &m) const
std::vector< MX > constraints() const
MX g() const
Get all (scalarised) constraint expressions as a column vector.
std::string describe(const MX &x, casadi_int indent=0) const
MetaCon get_meta_con(const MX &m) const
Get meta-data of symbol (for internal use only)
OptiAdvanced baked_copy() const
MX dual(const MX &m) const
get the dual variable
MX x() const
Get all (scalarised) decision variables as a symbolic column vector.
casadi_int nx() const
Number of (scalarised) decision variables.
MX variable(casadi_int n=1, casadi_int m=1, const std::string &attribute="full")
Create a decision variable (symbol)
void bake()
Fix the structure of the optimization problem.
casadi_int instance_number() const
DM value(const MX &x, const std::vector< MX > &values=std::vector< MX >()) const
void disp(std::ostream &stream, bool more=false) const override
Print representation.
std::vector< MX > active_symvar(VariableType type) const
void set_domain(const MX &x, const std::string &domain)
Set domain of variable.
static OptiNode * create(const std::string &problem_type)
OptiSol solve(bool accept_limit)
Crunch the numbers; solve the problem.
void res(const DMDict &res)
casadi_int np() const
Number of (scalarised) parameters.
bool has_callback_class() const
MX g_lookup(casadi_int i) const
void subject_to(const MX &g)
brief Add constraints
casadi_int ng() const
Number of (scalarised) constraints.
MetaVar get_meta(const MX &m) const
Get meta-data of symbol (for internal use only)
void subject_to()
Clear constraints.
void minimize(const MX &f)
Set objective.
void assert_baked() const
void update_user_dict(const MX &m, const Dict &meta)
add meta-data of an expression
void assert_solved() const
std::vector< MX > initial() const
get assignment expressions for initial values
void set_meta(const MX &m, const MetaVar &meta)
Set meta-data of an expression.
std::vector< DM > active_values(VariableType type) const
std::string class_name() const override
MX parameter(casadi_int n=1, casadi_int m=1, const std::string &attribute="full")
Create a parameter (symbol); fixed during optimization.
void set_initial(const std::vector< MX > &assignments)
DM value(const DM &x, const std::vector< MX > &values=std::vector< MX >()) const
bool is_parametric(const MX &expr) const
return true if expression is only dependant on Opti parameters, not variables
void assert_active_symbol(const MX &m) const
bool return_success(bool accept_limit) const
Did the solver return successfully?
bool problem_dirty() const
Function to_function(const std::string &name, const std::vector< MX > &args, const std::vector< MX > &res, const std::vector< std::string > &name_in, const std::vector< std::string > &name_out, const Dict &opts)
Create a CasADi Function from the Opti solver.
std::vector< MX > symvar() const
void callback_class(OptiCallback *callback)
std::vector< MX > value_variables() const
get assignment expressions for latest values
Opti copy() const
Copy.
MX x_lookup(casadi_int i) const
std::vector< MX > value_parameters() const
void mark_problem_dirty(bool flag=true)
OptiNode(const std::string &problem_type)
Create Opti Context.
std::string return_status() const
Get return status of solver.
~OptiNode()
Destructor.
std::vector< MX > symvar(const MX &expr, VariableType type) const
void assert_empty() const
void set_initial(const MX &x, const DM &v)
void mark_solver_dirty(bool flag=true)
MX p() const
Get all (scalarised) parameters as a symbolic column vector.
void set_value(const MX &x, const DM &v)
Set value of parameter.
MX lam_g() const
Get dual variables as a symbolic column vector.
MetaCon canon_expr(const MX &expr) const
Interpret an expression (for internal use only)
void set_value(const std::vector< MX > &assignments)
Set value of parameter.
DM value(const SX &x, const std::vector< MX > &values=std::vector< MX >()) const
std::string x_describe(casadi_int i) const
void solver(const std::string &solver, const Dict &plugin_options=Dict(), const Dict &solver_options=Dict())
Solver.
MX f() const
Get objective expression.
void mark_solved(bool flag=true)
Function casadi_solver() const
Get the underlying CasADi solver of the Opti stack.
std::vector< MX > symvar(const MX &expr) const
DMDict solve_actual(const DMDict &args)
Dict stats() const
Get statistics.
void set_meta_con(const MX &m, const MetaCon &meta)
Set meta-data of an expression.
A simplified interface for NLP modeling/solving.
Definition: optistack.hpp:628
A simplified interface for NLP modeling/solving.
Definition: optistack.hpp:90
The casadi namespace.
std::map< std::string, MX > MXDict
Definition: mx.hpp:948
GenericType::Dict Dict
C++ equivalent of Python's dict or MATLAB's struct.
std::map< std::string, DM > DMDict
Definition: dm_fwd.hpp:36