humoto
solver.h
Go to the documentation of this file.
1 /**
2  @file
3  @author Alexander Sherikov
4 
5  @brief
6 */
7 
8 #pragma once
9 
10 
11 namespace humoto
12 {
13  namespace qpmad
14  {
15  /**
16  * @brief Parameters of the solver
17  */
19  {
20  #define HUMOTO_CONFIG_SECTION_ID "SolverParameters"
21  #define HUMOTO_CONFIG_CONSTRUCTOR SolverParameters
22  #define HUMOTO_CONFIG_ENTRIES \
23  HUMOTO_CONFIG_PARENT_CLASS(SolverParametersBase) \
24  HUMOTO_CONFIG_SCALAR_(regularization_factor) \
25  HUMOTO_CONFIG_MEMBER_CLASS(options_, "QPmadParameters")
26  #include HUMOTO_CONFIG_DEFINE_ACCESSORS
27 
28 
29  public:
32  {
33  #define HUMOTO_CONFIG_SECTION_ID "QPmadParameters"
34  #define HUMOTO_CONFIG_CONSTRUCTOR QPmadParameters
35  #define HUMOTO_CONFIG_ENTRIES \
36  HUMOTO_CONFIG_ENUM_(hessian_type)\
37  HUMOTO_CONFIG_SCALAR_(tolerance)\
38  HUMOTO_CONFIG_SCALAR_(max_iter)
39  #include HUMOTO_CONFIG_DEFINE_ACCESSORS
40 
41  protected:
42  void setDefaults()
43  {
44  }
45 
46  public:
48  {
49  }
50  };
51 
52 
53  /// Regularization factor (disabled if = 0)
55 
57 
58 
59  public:
60  void setDefaults()
61  {
62  regularization_factor_ = 0;
63  }
64 
65  /**
66  * @brief Default constructor
67  */
69  {
70  setDefaults();
71  }
72  };
73 
74 
75 
76  /**
77  * @brief QP solver.
78  */
79  class HUMOTO_LOCAL Solver : public humoto::QPSolverMixin< humoto::Solver<SolverParameters> >
80  {
81  private:
83 
85 
86 
87  private:
88  /// @copydoc humoto::Solver::initialize
89  void initialize( const humoto::OptimizationProblem &hierarchy,
90  const humoto::SolutionStructure &sol_structure)
91  {
92  hierarchy.getQPProblem(qp_problem_, sol_structure, false);
93  }
94 
95 
96 
97  /// @copydoc humoto::Solver::solveHierarchy
99  const humoto::OptimizationProblem &hierarchy)
100  {
101  HUMOTO_ASSERT ( parameters_.regularization_factor_ >= 0,
102  "Regularization factor must be nonnegative.");
103 
104  if (parameters_.regularization_factor_ > 0)
105  {
106  qp_problem_.regularize(parameters_.regularization_factor_);
107  }
108 
109  ::qpmad::Solver::ReturnStatus qp_status =
110  solver.solve( solution.x_,
111  qp_problem_.getHessian(),
112  qp_problem_.getGradient(),
113  qp_problem_.getSolutionLowerBounds(),
114  qp_problem_.getSolutionUpperBounds(),
115  qp_problem_.getGeneralConstraints().getA(),
116  qp_problem_.getGeneralConstraints().getLowerBounds(),
117  qp_problem_.getGeneralConstraints().getUpperBounds());
118 
119 
120  switch (qp_status)
121  {
122  case ::qpmad::Solver::OK:
123  solution.return_status_ = SolverStatus::OK;
124  break;
125  case ::qpmad::Solver::MAXIMAL_NUMBER_OF_ITERATIONS:
127  break;
128  default:
130  break;
131  }
132 
133 
134  try
135  {
136  humoto::qpmad::Solution &qpmad_solution = dynamic_cast<humoto::qpmad::Solution &> (solution);
137 
138  qpmad_solution.qpmad_return_value_ = qp_status;
139  }
140  catch(...)
141  {
142  // not critical
143  }
144  }
145 
146 
147  public:
148  /**
149  * @brief Default constructor (with default parameters)
150  */
152  {
153  }
154 
155 
156  /**
157  * @brief Destructor
158  */
160  {
161  }
162 
163 
164  /**
165  * @brief Set parameters.
166  *
167  * @param[in] parameters parameters
168  */
169  explicit Solver(const SolverParameters &parameters)
170  {
171  setParameters(parameters);
172  }
173 
174 
175  /**
176  * @brief Log a QP problem
177  *
178  * @param[in,out] logger logger
179  * @param[in] parent parent
180  * @param[in] name name
181  */
183  const LogEntryName &parent = LogEntryName(),
184  const std::string &name = "qpmad") const
185  {
186  LogEntryName subname = parent; subname.add(name);
187  qp_problem_.log(logger, subname);
188  }
189  };
190  }
191 }
humoto::QPProblem_ILU_ALU qp_problem_
Definition: solver.h:84
const humoto::constraints::ContainerALU & getGeneralConstraints() const
get container with general constraints
Definition: qp_problem.h:250
~Solver()
Destructor.
Definition: solver.h:159
void solveHierarchy(humoto::Solution &solution, const humoto::OptimizationProblem &hierarchy)
Solve the hierarchy.
Definition: solver.h:98
#define HUMOTO_LOCAL
Definition: export_import.h:26
Eigen::VectorXd x_
Definition: solution.h:186
#define HUMOTO_GLOBAL_LOGGER_IF_DEFINED
Definition: logger.h:997
#define HUMOTO_ASSERT(condition, message)
void setDefaults()
Set members to their default values.
Definition: solver.h:42
Solution of a QP.
Definition: solution.h:19
Default configurable base is strict.
Definition: config.h:353
QPmadParameters options_
Definition: solver.h:56
Analog of &#39;sol_structure&#39; struct in Octave code. [determine_solution_structure.m].
Definition: solution.h:33
Container of the solution.
Definition: solution.h:176
Represents log entry name.
Definition: logger.h:169
Solver(const SolverParameters &parameters)
Set parameters.
Definition: solver.h:169
const Eigen::VectorXd & getSolutionLowerBounds() const
Get LB vector: &#39;LB <= X&#39;.
Definition: qp_problem.h:260
const Eigen::VectorXd & getGradient() const
Const accessor.
Definition: qp_problem.h:108
SolverStatus::Status return_status_
Definition: solution.h:187
Eigen::VectorXd & getUpperBounds()
Get upper bounds (ub vectors from &#39;A*x <= ub&#39;).
Eigen::MatrixXd & getA()
Get matrix A from general constraints: &#39;A*x = b&#39;, &#39;lb <= A*x <= ub&#39;.
An optimization problem [initialize_stack.m, simulation_loop.m].
Definition: hierarchy.h:20
::qpmad::Solver solver
Definition: solver.h:82
Threaded logger: any data sent to this logger is wrapped in a message and pushed to a queue...
Definition: logger.h:555
const Eigen::VectorXd & getSolutionUpperBounds() const
Get UB vector: &#39;X <= UB&#39;.
Definition: qp_problem.h:270
QP solver.
Definition: solver.h:79
A QP problem with constraints of a specific type.
Definition: qp_problem.h:450
void initialize(const humoto::OptimizationProblem &hierarchy, const humoto::SolutionStructure &sol_structure)
Initialize solver.
Definition: solver.h:89
The root namespace of HuMoTo.
Definition: config.h:12
void log(humoto::Logger &logger, const LogEntryName &parent=LogEntryName(), const std::string &name="qpmad") const
Log a QP problem.
Definition: solver.h:182
Parameters of the solver.
Definition: solver.h:18
Mixin solver interface (QP solver)
Definition: solver.h:462
Parameters of a solver.
Definition: solver.h:17
Solver()
Default constructor (with default parameters)
Definition: solver.h:151
void regularize(const double regularization_factor)
Adds regularization to the Hessian.
Definition: qp_problem.h:168
double regularization_factor_
Regularization factor (disabled if = 0)
Definition: solver.h:54
LogEntryName & add(const char *name)
extends entry name with a subname
Definition: logger.h:232
SolverParameters()
Default constructor.
Definition: solver.h:68
void solve(humoto::Solution &solution, const humoto::OptimizationProblem &hierarchy)
Solve an optimization problem (using previously specified parameters)
Definition: solver.h:155
Eigen::VectorXd & getLowerBounds()
Get lower bounds (lb/ub vectors from &#39;lb <= A*x <= ub&#39;).
void log(humoto::Logger &logger, const LogEntryName &parent=LogEntryName(), const std::string &name="qp_problem") const
Log a QP problem.
Definition: qp_problem.h:434
void setDefaults()
Set members to their default values.
Definition: solver.h:60
const Eigen::MatrixXd & getHessian() const
Const accessor.
Definition: qp_problem.h:98
void getQPProblem(QPProblemBase< t_QPConstraints > &qp_problem, const humoto::SolutionStructure &sol_structure, const bool initialize_upper_triangular_part=true) const
Form a QP problem.
Definition: hierarchy.h:329
::qpmad::Solver::ReturnStatus qpmad_return_value_
The return value of qpmad.
Definition: solution.h:23