optistack.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_HPP
26 #define CASADI_OPTISTACK_HPP
27 
28 #include "function.hpp"
29 #include "callback.hpp"
30 
31 namespace casadi {
32 
33 typedef DM native_DM;
34 
35 class OptiNode;
36 class OptiSol;
37 class OptiAdvanced;
38 class OptiCallback;
88 class CASADI_EXPORT Opti
89  : public SWIG_IF_ELSE(PrintableCommon, Printable<Opti>),
90  public SharedObject {
91  friend class InternalOptiCallback;
92 public:
93 
99  Opti(const std::string& problem_type="nlp");
100 
112  MX variable(casadi_int n=1, casadi_int m=1, const std::string& attribute="full");
113  MX variable(const Sparsity& sp, const std::string& attribute="full");
114  MX variable(const MX& symbol, const std::string& attribute="full");
115 
127  MX parameter(casadi_int n=1, casadi_int m=1, const std::string& attribute="full");
128  MX parameter(const Sparsity& sp, const std::string& attribute="full");
129  MX parameter(const MX& symbol, const std::string& attribute="full");
130 
137  void minimize(const MX& f, double linear_scale=1);
138 
140 
165  void subject_to(const MX& g, const Dict& options=Dict());
166  void subject_to(const std::vector<MX>& g, const Dict& options=Dict());
167  void subject_to(const MX& g, const DM& linear_scale, const Dict& options=Dict());
168  void subject_to(const std::vector<MX>& g, const DM& linear_scale, const Dict& options=Dict());
170 
172  void subject_to();
173 
184  void solver(const std::string& solver,
185  const Dict& plugin_options=Dict(),
186  const Dict& solver_options=Dict());
187 
189 
195  void set_initial(const MX& x, const DM& v);
196  void set_initial(const std::vector<MX>& assignments);
198 
200 
205  void set_value(const MX& x, const DM& v);
206  void set_value(const std::vector<MX>& assignments);
208 
210 
221  void set_domain(const MX& x, const std::string& domain);
222 
237  void set_linear_scale(const MX& x, const DM& scale, const DM& offset=0);
238 
241 
249 
251 
259  native_DM value(const MX& x, const std::vector<MX>& values=std::vector<MX>()) const;
260  native_DM value(const DM& x, const std::vector<MX>& values=std::vector<MX>()) const;
261  native_DM value(const SX& x, const std::vector<MX>& values=std::vector<MX>()) const;
263 
270  Dict stats() const;
271 
278  std::string return_status() const;
279 
283  std::vector<MX> initial() const;
284 
288  std::vector<MX> value_variables() const;
289  std::vector<MX> value_parameters() const;
290 
294  Function scale_helper(const Function& h) const;
295 
303  MX dual(const MX& m) const;
304 
308  casadi_int nx() const;
309 
313  casadi_int np() const;
314 
318  casadi_int ng() const;
319 
323  MX x() const;
324 
328  MX p() const;
329 
333  MX g() const;
334 
338  MX f() const;
339 
343  MX lbg() const;
344  MX ubg() const;
345 
349  double f_linear_scale() const;
350 
360  MX lam_g() const;
361 
363 
372  Function to_function(const std::string& name,
373  const std::vector<MX>& args, const std::vector<MX>& res,
374  const Dict& opts = Dict());
375 
376  Function to_function(const std::string& name,
377  const std::vector<MX>& args, const std::vector<MX>& res,
378  const std::vector<std::string>& name_in,
379  const std::vector<std::string>& name_out,
380  const Dict& opts = Dict());
381 
382  Function to_function(const std::string& name,
383  const std::map<std::string, MX>& dict,
384  const std::vector<std::string>& name_in,
385  const std::vector<std::string>& name_out,
386  const Dict& opts = Dict());
388 
389  #ifndef SWIGMATLAB
397  static MX bounded(const MX& lb, const MX& expr, const MX& ub) { return (lb<=expr)<= ub; }
398  #endif
399 
409 
419 
426  Opti copy() const;
427 
434  void update_user_dict(const MX& m, const Dict& meta);
435  void update_user_dict(const std::vector<MX>& m, const Dict& meta);
437  Dict user_dict(const MX& m) const;
438 
440  std::string type_name() const { return "Opti"; }
441 
443  void disp(std::ostream& stream, bool more=false) const;
444 
446  std::string get_str(bool more=false) const;
447 
449 
454  void callback_class(OptiCallback* callback);
457 
458 #ifndef SWIG
459  Opti(const Opti& x);
460 
464  ~Opti() {}
465 
466  static Opti create(OptiNode* node);
469 
472  OptiNode* operator->();
473 
477  const OptiNode* operator->() const;
480 
481  Opti(OptiNode* node);
482 
483 #endif // SWIG
484 
485 };
486 
488  OPTI_GENERIC_EQUALITY, // g1(x,p) == g2(x,p)
489  OPTI_GENERIC_INEQUALITY, // g1(x,p) <= g2(x,p)
490  OPTI_EQUALITY, // g(x,p) == bound(p)
491  OPTI_INEQUALITY, // g(x,p) <= bound(p)
492  OPTI_DOUBLE_INEQUALITY, // lb(p) <= g(x,p) <= ub(p)
493  OPTI_PSD, // A(x,p) >= b(p)
496  OPTI_VAR, // variable
497  OPTI_PAR, // parameter
498  OPTI_DUAL_G // dual
499  };
500  enum DomainType {
503  };
504 
505 
508  casadi_int start;
509  casadi_int stop;
510  };
512  MetaCon() : n(1), flipped(false), linear_scale(1) {}
513  MX original; // original expression
514  MX canon; // Canonical expression
518  casadi_int n;
519  bool flipped;
524  };
526  std::string attribute;
527  casadi_int n;
528  casadi_int m;
531  casadi_int count;
532  casadi_int i;
533  casadi_int active_i;
535  };
536 
537 
539 public:
541  }
543  casadi_error("Callback objects cannot be copied");
544  }
545  virtual void call(casadi_int i) {
546  uout() << "This is a simple callback at iteration" << i << std::endl;
547  }
548  virtual ~OptiCallback() {}
549 };
550 
551 class CASADI_EXPORT OptiAdvanced : public Opti {
552  friend class InternalOptiCallback;
553 public:
554 
555  OptiAdvanced(const Opti& x);
556 
561 
562 
565 
567  bool is_parametric(const MX& expr) const;
568 
570 
576  std::vector<MX> symvar() const;
577  std::vector<MX> symvar(const MX& expr) const;
578  std::vector<MX> symvar(const MX& expr, VariableType type) const;
580 
582  MetaCon canon_expr(const MX& expr) const;
583 
585  MetaVar get_meta(const MX& m) const;
586 
588  MetaCon get_meta_con(const MX& m) const;
589 
591  void set_meta(const MX& m, const MetaVar& meta);
592 
594  void set_meta_con(const MX& m, const MetaCon& meta);
595 
596  void assert_active_symbol(const MX& m) const;
597 
598  std::vector<MX> active_symvar(VariableType type) const;
599  std::vector<DM> active_values(VariableType type) const;
600 
601  MX x_lookup(casadi_index i) const;
602  MX g_lookup(casadi_index i) const;
603 
604  casadi_index g_index_reduce_g(casadi_index i) const;
605  casadi_index g_index_reduce_x(casadi_index i) const;
606  casadi_index g_index_unreduce_g(casadi_index i) const;
607 
608  std::string x_describe(casadi_index i, const Dict& opts=Dict()) const;
609  std::string g_describe(casadi_index i, const Dict& opts=Dict()) const;
610  std::string describe(const MX& x, casadi_index indent=0, const Dict& opts=Dict()) const;
611 
612  void show_infeasibilities(double tol=0, const Dict& opts=Dict()) const;
613 
615  DMDict solve_actual(const DMDict& args);
616 
617  DMDict arg() const;
618  void res(const DMDict& res);
619  DMDict res() const;
620  std::vector<MX> constraints() const;
621  MX objective() const;
622 
624 
625  void assert_empty() const;
626 
627 
629  void bake();
630 
632  void mark_problem_dirty(bool flag=true);
633  bool problem_dirty() const;
634 
636  void mark_solver_dirty(bool flag=true);
637  bool solver_dirty() const;
638 
639  bool solved_;
640  void mark_solved(bool flag=true);
641  bool solved() const;
642 
643  void assert_solved() const;
644  void assert_baked() const;
645 
646  casadi_int instance_number() const;
647 
648 protected:
650 };
651 
662 class CASADI_EXPORT OptiSol : public SWIG_IF_ELSE(PrintableCommon, Printable<OptiAdvanced>) {
663  friend class OptiNode;
664  public:
665  std::string type_name() const {return "OptiSol";}
666  void disp(std::ostream& stream, bool more=false) const;
667  std::string get_str(bool more=false) const;
669 
677  native_DM value(const MX& x, const std::vector<MX>& values=std::vector<MX>()) const;
678  native_DM value(const DM& x, const std::vector<MX>& values=std::vector<MX>()) const;
679  native_DM value(const SX& x, const std::vector<MX>& values=std::vector<MX>()) const;
681 
683  std::vector<MX> value_variables() const;
684  std::vector<MX> value_parameters() const;
685 
692  Dict stats() const;
693 
694  Opti opti() const { return optistack_; } // NOLINT(cppcoreguidelines-slicing)
695 
696  protected:
697  OptiSol(const Opti& opti);
699 };
700 
701 
702 } // namespace casadi
703 
704 #endif // CASADI_OPTI_HPP
Function object.
Definition: function.hpp:60
MX - Matrix expression.
Definition: mx.hpp:92
casadi_index g_index_unreduce_g(casadi_index i) const
void mark_problem_dirty(bool flag=true)
void mark_solver_dirty(bool flag=true)
MX x_lookup(casadi_index i) const
bool problem_dirty() const
std::string describe(const MX &x, casadi_index indent=0, const Dict &opts=Dict()) const
std::vector< MX > symvar(const MX &expr) const
Get symbols present in expression.
std::vector< MX > constraints() const
DMDict solve_actual(const DMDict &args)
void bake()
Fix the structure of the optimization problem.
std::vector< MX > active_symvar(VariableType type) const
MX g_lookup(casadi_index i) const
OptiAdvanced(const Opti &x)
bool solver_dirty() const
void set_meta_con(const MX &m, const MetaCon &meta)
Set meta-data of an expression.
casadi_index g_index_reduce_x(casadi_index i) const
MetaCon canon_expr(const MX &expr) const
Interpret an expression (for internal use only)
std::string x_describe(casadi_index i, const Dict &opts=Dict()) const
std::vector< MX > symvar() const
Get symbols present in expression.
Function casadi_solver() const
Get the underlying CasADi solver of the Opti stack.
void res(const DMDict &res)
std::vector< MX > symvar(const MX &expr, VariableType type) const
Get symbols present in expression.
void assert_baked() const
std::string g_describe(casadi_index i, const Dict &opts=Dict()) const
void assert_empty() const
casadi_index g_index_reduce_g(casadi_index i) const
void mark_solved(bool flag=true)
OptiAdvanced baked_copy() const
bool is_parametric(const MX &expr) const
return true if expression is only dependant on Opti parameters, not variables
DMDict res() const
void show_infeasibilities(double tol=0, const Dict &opts=Dict()) const
std::vector< DM > active_values(VariableType type) const
void assert_active_symbol(const MX &m) const
void assert_solved() const
DMDict arg() const
bool solved() const
~OptiAdvanced()
Destructor.
Definition: optistack.hpp:560
void set_meta(const MX &m, const MetaVar &meta)
Set meta-data of an expression.
MetaVar get_meta(const MX &m) const
Get meta-data of symbol (for internal use only)
MetaCon get_meta_con(const MX &m) const
Get meta-data of symbol (for internal use only)
casadi_int instance_number() const
virtual ~OptiCallback()
Definition: optistack.hpp:548
virtual void call(casadi_int i)
Definition: optistack.hpp:545
OptiCallback(const OptiCallback &obj)
Definition: optistack.hpp:542
A simplified interface for NLP modeling/solving.
A simplified interface for NLP modeling/solving.
Definition: optistack.hpp:662
std::string get_str(bool more=false) const
std::vector< MX > value_parameters() const
void disp(std::ostream &stream, bool more=false) const
native_DM value(const DM &x, const std::vector< MX > &values=std::vector< MX >()) const
std::string type_name() const
Definition: optistack.hpp:665
Dict stats() const
Get statistics.
Opti opti() const
Definition: optistack.hpp:694
native_DM value(const MX &x, const std::vector< MX > &values=std::vector< MX >()) const
native_DM value(const SX &x, const std::vector< MX > &values=std::vector< MX >()) const
OptiAdvanced optistack_
Definition: optistack.hpp:698
std::vector< MX > value_variables() const
get assignment expressions for the optimal solution
OptiSol(const Opti &opti)
A simplified interface for NLP modeling/solving.
Definition: optistack.hpp:90
native_DM value(const MX &x, const std::vector< MX > &values=std::vector< MX >()) const
native_DM value(const SX &x, const std::vector< MX > &values=std::vector< MX >()) const
Set domain of a decision variable.
std::vector< MX > value_variables() const
get assignment expressions for latest values
Function scale_helper(const Function &h) const
Scale a helper function constructed via opti.x, opti.g, ...
void set_domain(const MX &x, const std::string &domain)
Set domain of a decision variable.
void subject_to(const MX &g, const DM &linear_scale, const Dict &options=Dict())
Add constraints.
DM x_linear_scale() const
Function to_function(const std::string &name, const std::map< std::string, MX > &dict, const std::vector< std::string > &name_in, const std::vector< std::string > &name_out, const Dict &opts=Dict())
Create a CasADi Function from the Opti solver.
native_DM value(const DM &x, const std::vector< MX > &values=std::vector< MX >()) const
Set domain of a decision variable.
void subject_to()
Clear constraints.
MX parameter(const Sparsity &sp, const std::string &attribute="full")
static MX bounded(const MX &lb, const MX &expr, const MX &ub)
Construct a double inequality.
Definition: optistack.hpp:397
casadi_int ng() const
Number of (scalarised) constraints.
void minimize(const MX &f, double linear_scale=1)
Set objective.
void callback_class(OptiCallback *callback)
Helper methods for callback()
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=Dict())
Create a CasADi Function from the Opti solver.
casadi_int np() const
Number of (scalarised) parameters.
void set_initial(const std::vector< MX > &assignments)
Dict stats() const
Get statistics.
void subject_to(const std::vector< MX > &g, const Dict &options=Dict())
Add constraints.
MX variable(const MX &symbol, const std::string &attribute="full")
MX parameter(const MX &symbol, const std::string &attribute="full")
void subject_to(const MX &g, const Dict &options=Dict())
Add constraints.
void disp(std::ostream &stream, bool more=false) const
Print representation.
MX lbg() const
Get all (scalarised) bounds on constraints as a column vector.
std::string return_status() const
Get return status of solver.
casadi_int nx() const
Number of (scalarised) decision variables.
Opti(const std::string &problem_type="nlp")
Create Opti Context.
std::string get_str(bool more=false) const
Get string representation.
OptiAdvanced debug() const
Get a copy with advanced functionality.
MX variable(casadi_int n=1, casadi_int m=1, const std::string &attribute="full")
Create a decision variable (symbol)
MX f() const
Get objective expression.
OptiAdvanced advanced() const
Get a copy with advanced functionality.
Opti copy() const
Get a copy of the.
Dict user_dict(const MX &m) const
Get user data.
void subject_to(const std::vector< MX > &g, const DM &linear_scale, const Dict &options=Dict())
Add constraints.
Function to_function(const std::string &name, const std::vector< MX > &args, const std::vector< MX > &res, const Dict &opts=Dict())
Create a CasADi Function from the Opti solver.
void set_linear_scale(const MX &x, const DM &scale, const DM &offset=0)
Set scale of a decision variable.
std::vector< MX > initial() const
get assignment expressions for initial values
MX parameter(casadi_int n=1, casadi_int m=1, const std::string &attribute="full")
Create a parameter (symbol); fixed during optimization.
void set_value(const MX &x, const DM &v)
Set value of parameter.
DM x_linear_scale_offset() const
void callback_class()
Helper methods for callback()
MX p() const
Get all (scalarised) parameters as a symbolic column vector.
MX variable(const Sparsity &sp, const std::string &attribute="full")
MX ubg() const
void solver(const std::string &solver, const Dict &plugin_options=Dict(), const Dict &solver_options=Dict())
Set a solver.
std::vector< MX > value_parameters() const
MX dual(const MX &m) const
get the dual variable
void set_initial(const MX &x, const DM &v)
void update_user_dict(const MX &m, const Dict &meta)
add user data
OptiSol solve()
Crunch the numbers; solve the problem.
void set_value(const std::vector< MX > &assignments)
Set value of parameter.
DM g_linear_scale() const
MX lam_g() const
Get all (scalarised) dual variables as a symbolic column vector.
MX x() const
Get all (scalarised) decision variables as a symbolic column vector.
void update_user_dict(const std::vector< MX > &m, const Dict &meta)
std::string type_name() const
Readable name of the class.
Definition: optistack.hpp:440
double f_linear_scale() const
OptiSol solve_limited()
Crunch the numbers; solve the problem.
MX g() const
Get all (scalarised) constraint expressions as a column vector.
GenericShared implements a reference counting framework similar for efficient and.
General sparsity class.
Definition: sparsity.hpp:106
The casadi namespace.
Definition: archiver.hpp:32
CASADI_EXPORT std::ostream & uout()
@ OPTI_VAR
Definition: optistack.hpp:496
@ OPTI_DUAL_G
Definition: optistack.hpp:498
@ OPTI_PAR
Definition: optistack.hpp:497
@ OPTI_DOMAIN_INTEGER
Definition: optistack.hpp:502
@ OPTI_DOMAIN_REAL
Definition: optistack.hpp:501
GenericType::Dict Dict
C++ equivalent of Python's dict or MATLAB's struct.
DM native_DM
Definition: optistack.hpp:33
ConstraintType
Definition: optistack.hpp:487
@ OPTI_DOUBLE_INEQUALITY
Definition: optistack.hpp:492
@ OPTI_INEQUALITY
Definition: optistack.hpp:491
@ OPTI_GENERIC_INEQUALITY
Definition: optistack.hpp:489
@ OPTI_EQUALITY
Definition: optistack.hpp:490
@ OPTI_GENERIC_EQUALITY
Definition: optistack.hpp:488
@ OPTI_UNKNOWN
Definition: optistack.hpp:494
@ OPTI_PSD
Definition: optistack.hpp:493
std::map< std::string, DM > DMDict
Definition: dm_fwd.hpp:36
ConstraintType type
Definition: optistack.hpp:515
casadi_int n
Definition: optistack.hpp:518
casadi_int active_i
Definition: optistack.hpp:533
std::string attribute
Definition: optistack.hpp:526
DomainType domain
Definition: optistack.hpp:530
VariableType type
Definition: optistack.hpp:529
casadi_int count
Definition: optistack.hpp:531
casadi_int i
Definition: optistack.hpp:532
casadi_int m
Definition: optistack.hpp:528
casadi_int n
Definition: optistack.hpp:527