robotoc
robotoc - efficient ROBOT Optimal Control solvers
Loading...
Searching...
No Matches
sto_constraints.hpp
Go to the documentation of this file.
1#ifndef ROBOTOC_STO_CONSTRAINTS_HPP_
2#define ROBOTOC_STO_CONSTRAINTS_HPP_
3
4#include <limits>
5#include <cmath>
6#include <vector>
7
8#include "Eigen/Core"
9
15
16
17namespace robotoc {
18
24public:
37 STOConstraints(const int num_switches,
38 const double minimum_dwell_time=std::sqrt(std::numeric_limits<double>::epsilon()),
39 const double barrier_param=1.0e-03,
40 const double fraction_to_boundary_rule=0.995);
41
52 STOConstraints(const std::vector<double>& minimum_dwell_times,
53 const double barrier_param=1.0e-03,
54 const double fraction_to_boundary_rule=0.995);
55
66 STOConstraints(const Eigen::VectorXd& minimum_dwell_times,
67 const double barrier_param=1.0e-03,
68 const double fraction_to_boundary_rule=0.995);
69
74
78 ~STOConstraints() = default;
79
83 STOConstraints(const STOConstraints&) = default;
84
89
93 STOConstraints(STOConstraints&&) noexcept = default;
94
98 STOConstraints& operator=(STOConstraints&&) noexcept = default;
99
107 const double minimum_dwell_time=std::sqrt(std::numeric_limits<double>::epsilon()));
108
114 void setMinimumDwellTimes(const std::vector<double>& minimum_dwell_times);
115
121 void setMinimumDwellTimes(const Eigen::VectorXd& minimum_dwell_times);
122
127 const Eigen::VectorXd& getMinimumDwellTimes() const;
128
133 void setBarrierParam(const double barrier_param);
134
141 void setFractionToBoundaryRule(const double fraction_to_boundary_rule);
142
147 double getBarrierParam() const;
148
154
161 const TimeDiscretization& time_discretization) const;
162
169 bool isFeasible(const TimeDiscretization& time_discretization,
170 ConstraintComponentData& data) const;
171
177 void setSlackAndDual(const TimeDiscretization& time_discretization,
178 ConstraintComponentData& data) const;
179
186 void evalConstraint(const TimeDiscretization& time_discretization,
187 ConstraintComponentData& data) const;
188
197 void linearizeConstraints(const TimeDiscretization& time_discretization,
199 Eigen::VectorXd& lt) const;
200
210 void condenseSlackAndDual(ConstraintComponentData& data, Eigen::VectorXd& lt,
211 Eigen::MatrixXd& Qtt) const;
212
219 void expandSlackAndDual(ConstraintComponentData& data, Eigen::VectorXd& dts) const;
220
227 double maxSlackStepSize(const ConstraintComponentData& data) const;
228
235 double maxDualStepSize(const ConstraintComponentData& data) const;
236
242 void updateSlack(ConstraintComponentData& data, const double step_size) const;
243
249 void updateDual(ConstraintComponentData& data, const double step_size) const;
250
256 static void computeDwellTimes(const TimeDiscretization& time_discretization,
257 Eigen::VectorXd& dwell_times);
258
259private:
260 double barrier_, fraction_to_boundary_rule_;
261 int num_switches_;
262
263 Eigen::VectorXd primal_step_size_, dual_step_size_;
264 Eigen::VectorXd minimum_dwell_times_;
265
266};
267
268} // namespace robotoc
269
270#endif // ROBOTOC_STO_CONSTRAINTS_HPP_
Data used in constraint components. Composed by slack, dual (Lagrange multiplier),...
Definition: constraint_component_data.hpp:17
Constraints of the switching time optimization problems.
Definition: sto_constraints.hpp:23
void setSlackAndDual(const TimeDiscretization &time_discretization, ConstraintComponentData &data) const
Sets the slack and dual variables of constraints.
const Eigen::VectorXd & getMinimumDwellTimes() const
Gets the minimum dwell times.
STOConstraints(const STOConstraints &)=default
Default copy constructor.
void linearizeConstraints(const TimeDiscretization &time_discretization, ConstraintComponentData &data, Eigen::VectorXd &lt) const
Evaluates the constraints (i.e., calls evalConstraint()) and adds the products of the Jacobian of the...
STOConstraints(const int num_switches, const double minimum_dwell_time=std::sqrt(std::numeric_limits< double >::epsilon()), const double barrier_param=1.0e-03, const double fraction_to_boundary_rule=0.995)
Constructor.
double maxDualStepSize(const ConstraintComponentData &data) const
Computes and returns the maximum step size by applying fraction-to-boundary-rule to the directions of...
void evalConstraint(const TimeDiscretization &time_discretization, ConstraintComponentData &data) const
Computes the primal residual, residual in the complementary slackness, and the log-barrier function o...
void updateDual(ConstraintComponentData &data, const double step_size) const
Updates the dual variables according to step_size.
STOConstraints(const Eigen::VectorXd &minimum_dwell_times, const double barrier_param=1.0e-03, const double fraction_to_boundary_rule=0.995)
Constructor.
double getBarrierParam() const
Gets the barrier parameter.
double maxSlackStepSize(const ConstraintComponentData &data) const
Computes and returns the maximum step size by applying fraction-to-boundary-rule to the directions of...
void updateSlack(ConstraintComponentData &data, const double step_size) const
Updates the slack variables according to step_size.
void setMinimumDwellTimes(const double minimum_dwell_time=std::sqrt(std::numeric_limits< double >::epsilon()))
Sets the minimum dwell times.
static void computeDwellTimes(const TimeDiscretization &time_discretization, Eigen::VectorXd &dwell_times)
Computes the dwell times from the given time discretization.
bool isFeasible(const TimeDiscretization &time_discretization, ConstraintComponentData &data) const
Checks whether the current split solution s is feasible or not.
STOConstraints(STOConstraints &&) noexcept=default
Default move constructor.
STOConstraints(const std::vector< double > &minimum_dwell_times, const double barrier_param=1.0e-03, const double fraction_to_boundary_rule=0.995)
Constructor.
STOConstraints()
Default constructor.
ConstraintComponentData createConstraintsData(const TimeDiscretization &time_discretization) const
Creates the data for this constraints.
void setFractionToBoundaryRule(const double fraction_to_boundary_rule)
Sets the parameter of the fraction-to-boundary-rule for all the constraint components.
~STOConstraints()=default
Default destructor.
void condenseSlackAndDual(ConstraintComponentData &data, Eigen::VectorXd &lt, Eigen::MatrixXd &Qtt) const
Condenses the slack and dual variables. linearizeConstraints() must be called before this function.
void expandSlackAndDual(ConstraintComponentData &data, Eigen::VectorXd &dts) const
Expands the slack and dual, i.e., computes the directions of the slack and dual variables from the di...
void setBarrierParam(const double barrier_param)
Sets the barrier parameter for all the constraint components.
STOConstraints & operator=(const STOConstraints &)=default
Default copy operator.
double getFractionToBoundaryRule() const
Gets the parameter of the fraction-to-boundary-rule.
Time discretization of the optimal control problem.
Definition: time_discretization.hpp:20
Definition: constraint_component_base.hpp:17