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 
125  MX parameter(casadi_int n=1, casadi_int m=1, const std::string& attribute="full");
126 
133  void minimize(const MX& f);
134 
136 
161  void subject_to(const MX& g);
162  void subject_to(const std::vector<MX>& g);
164 
166  void subject_to();
167 
178  void solver(const std::string& solver,
179  const Dict& plugin_options=Dict(),
180  const Dict& solver_options=Dict());
181 
183 
189  void set_initial(const MX& x, const DM& v);
190  void set_initial(const std::vector<MX>& assignments);
192 
194 
199  void set_value(const MX& x, const DM& v);
200  void set_value(const std::vector<MX>& assignments);
202 
204 
215  void set_domain(const MX& x, const std::string& domain);
216 
218 
219 
222 
230 
232 
240  native_DM value(const MX& x, const std::vector<MX>& values=std::vector<MX>()) const;
241  native_DM value(const DM& x, const std::vector<MX>& values=std::vector<MX>()) const;
242  native_DM value(const SX& x, const std::vector<MX>& values=std::vector<MX>()) const;
244 
251  Dict stats() const;
252 
259  std::string return_status() const;
260 
264  std::vector<MX> initial() const;
265 
269  std::vector<MX> value_variables() const;
270  std::vector<MX> value_parameters() const;
271 
279  MX dual(const MX& m) const;
280 
284  casadi_int nx() const;
285 
289  casadi_int np() const;
290 
294  casadi_int ng() const;
295 
299  MX x() const;
300 
304  MX p() const;
305 
309  MX g() const;
310 
314  MX f() const;
315 
319  MX lbg() const;
320  MX ubg() const;
321 
331  MX lam_g() const;
332 
334 
343  Function to_function(const std::string& name,
344  const std::vector<MX>& args, const std::vector<MX>& res,
345  const Dict& opts = Dict());
346 
347  Function to_function(const std::string& name,
348  const std::vector<MX>& args, const std::vector<MX>& res,
349  const std::vector<std::string>& name_in,
350  const std::vector<std::string>& name_out,
351  const Dict& opts = Dict());
352 
353  Function to_function(const std::string& name,
354  const std::map<std::string, MX>& dict,
355  const std::vector<std::string>& name_in,
356  const std::vector<std::string>& name_out,
357  const Dict& opts = Dict());
359 
360  #ifndef SWIGMATLAB
368  static MX bounded(const MX& lb, const MX& expr, const MX& ub) { return (lb<=expr)<= ub; }
369  #endif
370 
380 
390 
397  Opti copy() const;
398 
405  void update_user_dict(const MX& m, const Dict& meta);
406  void update_user_dict(const std::vector<MX>& m, const Dict& meta);
408  Dict user_dict(const MX& m) const;
409 
411  std::string type_name() const { return "Opti"; }
412 
414  void disp(std::ostream& stream, bool more=false) const;
415 
417  std::string get_str(bool more=false) const;
418 
420 
425  void callback_class(OptiCallback* callback);
428 
429 #ifndef SWIG
430  Opti(const Opti& x);
431 
435  ~Opti() {}
436 
437  static Opti create(OptiNode* node);
440 
443  OptiNode* operator->();
444 
448  const OptiNode* operator->() const;
451 
452  Opti(OptiNode* node);
453 
454 #endif // SWIG
455 
456 };
457 
459  OPTI_GENERIC_EQUALITY, // g1(x,p) == g2(x,p)
460  OPTI_GENERIC_INEQUALITY, // g1(x,p) <= g2(x,p)
461  OPTI_EQUALITY, // g(x,p) == bound(p)
462  OPTI_INEQUALITY, // g(x,p) <= bound(p)
463  OPTI_DOUBLE_INEQUALITY, // lb(p) <= g(x,p) <= ub(p)
464  OPTI_PSD, // A(x,p) >= b(p)
467  OPTI_VAR, // variable
468  OPTI_PAR, // parameter
469  OPTI_DUAL_G // dual
470  };
471  enum DomainType {
474  };
475 
476 
479  casadi_int start;
480  casadi_int stop;
481  };
483  MetaCon() : n(1), flipped(false) {}
484  MX original; // original expression
485  MX canon; // Canonical expression
489  casadi_int n;
490  bool flipped;
494  };
496  std::string attribute;
497  casadi_int n;
498  casadi_int m;
501  casadi_int count;
502  casadi_int i;
503  casadi_int active_i;
505  };
506 
507 
509 public:
511  }
513  casadi_error("Callback objects cannot be copied");
514  }
515  virtual void call(casadi_int i) {
516  uout() << "This is a simple callback at iteration" << i << std::endl;
517  }
518  virtual ~OptiCallback() {}
519 };
520 
521 class CASADI_EXPORT OptiAdvanced : public Opti {
522  friend class InternalOptiCallback;
523 public:
524 
525  OptiAdvanced(const Opti& x);
526 
531 
532 
535 
537  bool is_parametric(const MX& expr) const;
538 
540 
546  std::vector<MX> symvar() const;
547  std::vector<MX> symvar(const MX& expr) const;
548  std::vector<MX> symvar(const MX& expr, VariableType type) const;
550 
552  MetaCon canon_expr(const MX& expr) const;
553 
555  MetaVar get_meta(const MX& m) const;
556 
558  MetaCon get_meta_con(const MX& m) const;
559 
561  void set_meta(const MX& m, const MetaVar& meta);
562 
564  void set_meta_con(const MX& m, const MetaCon& meta);
565 
566  void assert_active_symbol(const MX& m) const;
567 
568  std::vector<MX> active_symvar(VariableType type) const;
569  std::vector<DM> active_values(VariableType type) const;
570 
571  MX x_lookup(casadi_index i) const;
572  MX g_lookup(casadi_index i) const;
573 
574  std::string x_describe(casadi_index i) const;
575  std::string g_describe(casadi_index i) const;
576  std::string describe(const MX& x, casadi_index indent=0) const;
577 
578  void show_infeasibilities(double tol=0) const;
579 
581  DMDict solve_actual(const DMDict& args);
582 
583  DMDict arg() const;
584  void res(const DMDict& res);
585  DMDict res() const;
586  std::vector<MX> constraints() const;
587  MX objective() const;
588 
590 
591  void assert_empty() const;
592 
593 
595  void bake();
596 
598  void mark_problem_dirty(bool flag=true);
599  bool problem_dirty() const;
600 
602  void mark_solver_dirty(bool flag=true);
603  bool solver_dirty() const;
604 
605  bool solved_;
606  void mark_solved(bool flag=true);
607  bool solved() const;
608 
609  void assert_solved() const;
610  void assert_baked() const;
611 
612  casadi_int instance_number() const;
613 
614 protected:
616 };
617 
628 class CASADI_EXPORT OptiSol : public SWIG_IF_ELSE(PrintableCommon, Printable<OptiAdvanced>) {
629  friend class OptiNode;
630  public:
631  std::string type_name() const {return "OptiSol";}
632  void disp(std::ostream& stream, bool more=false) const;
633  std::string get_str(bool more=false) const;
635 
643  native_DM value(const MX& x, const std::vector<MX>& values=std::vector<MX>()) const;
644  native_DM value(const DM& x, const std::vector<MX>& values=std::vector<MX>()) const;
645  native_DM value(const SX& x, const std::vector<MX>& values=std::vector<MX>()) const;
647 
649  std::vector<MX> value_variables() const;
650  std::vector<MX> value_parameters() const;
651 
658  Dict stats() const;
659 
660  Opti opti() const { return optistack_; } // NOLINT(cppcoreguidelines-slicing)
661 
662  protected:
663  OptiSol(const Opti& opti);
665 };
666 
667 
668 } // namespace casadi
669 
670 #endif // CASADI_OPTI_HPP
Function object.
Definition: function.hpp:60
MX - Matrix expression.
Definition: mx.hpp:84
Sparse matrix class. SX and DM are specializations.
Definition: matrix_decl.hpp:92
void mark_problem_dirty(bool flag=true)
void mark_solver_dirty(bool flag=true)
void show_infeasibilities(double tol=0) const
MX x_lookup(casadi_index i) const
bool problem_dirty() const
std::string g_describe(casadi_index i) 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.
MetaCon canon_expr(const MX &expr) const
Interpret an expression (for internal use only)
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
void assert_empty() 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
std::string x_describe(casadi_index i) const
DMDict res() 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:530
std::string describe(const MX &x, casadi_index indent=0) const
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:518
virtual void call(casadi_int i)
Definition: optistack.hpp:515
OptiCallback(const OptiCallback &obj)
Definition: optistack.hpp:512
A simplified interface for NLP modeling/solving.
A simplified interface for NLP modeling/solving.
Definition: optistack.hpp:628
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:631
Dict stats() const
Get statistics.
Opti opti() const
Definition: optistack.hpp:660
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:664
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
std::vector< MX > value_variables() const
get assignment expressions for latest values
void set_domain(const MX &x, const std::string &domain)
Set domain of a decision variable.
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
void subject_to()
Clear constraints.
static MX bounded(const MX &lb, const MX &expr, const MX &ub)
Construct a double inequality.
Definition: optistack.hpp:368
casadi_int ng() const
Number of (scalarised) constraints.
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 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.
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.
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.
void callback_class()
Helper methods for callback()
void subject_to(const MX &g)
Add constraints.
MX p() const
Get all (scalarised) parameters as a symbolic column vector.
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 minimize(const MX &f)
Set objective.
void set_value(const std::vector< MX > &assignments)
Set value of parameter.
void subject_to(const std::vector< MX > &g)
Add constraints.
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:411
OptiSol solve_limited()
Crunch the numbers; solve the problem.
MX g() const
Get all (scalarised) constraint expressions as a column vector.
SharedObject implements a reference counting framework similar for efficient and.
The casadi namespace.
CASADI_EXPORT std::ostream & uout()
@ OPTI_VAR
Definition: optistack.hpp:467
@ OPTI_DUAL_G
Definition: optistack.hpp:469
@ OPTI_PAR
Definition: optistack.hpp:468
@ OPTI_DOMAIN_INTEGER
Definition: optistack.hpp:473
@ OPTI_DOMAIN_REAL
Definition: optistack.hpp:472
GenericType::Dict Dict
C++ equivalent of Python's dict or MATLAB's struct.
DM native_DM
Definition: optistack.hpp:33
ConstraintType
Definition: optistack.hpp:458
@ OPTI_DOUBLE_INEQUALITY
Definition: optistack.hpp:463
@ OPTI_INEQUALITY
Definition: optistack.hpp:462
@ OPTI_GENERIC_INEQUALITY
Definition: optistack.hpp:460
@ OPTI_EQUALITY
Definition: optistack.hpp:461
@ OPTI_GENERIC_EQUALITY
Definition: optistack.hpp:459
@ OPTI_UNKNOWN
Definition: optistack.hpp:465
@ OPTI_PSD
Definition: optistack.hpp:464
std::map< std::string, DM > DMDict
Definition: dm_fwd.hpp:36
ConstraintType type
Definition: optistack.hpp:486
casadi_int n
Definition: optistack.hpp:489
casadi_int active_i
Definition: optistack.hpp:503
std::string attribute
Definition: optistack.hpp:496
DomainType domain
Definition: optistack.hpp:500
VariableType type
Definition: optistack.hpp:499
casadi_int count
Definition: optistack.hpp:501
casadi_int i
Definition: optistack.hpp:502
casadi_int m
Definition: optistack.hpp:498
casadi_int n
Definition: optistack.hpp:497